#include #include #include #include "FixVec.h" #include "Boundary.h" #include "DOFAdmin.h" #include "ElInfo.h" #include "Error.h" #include "FiniteElemSpace.h" #include "Global.h" #include "Mesh.h" #include "Quadrature.h" #include "AbstractFunction.h" #include "BoundaryManager.h" #include "ElementVector.h" #include "Assembler.h" #include "OpenMP.h" #include "Operator.h" #include "Parameters.h" #include "Traverse.h" namespace AMDiS { template DOFVectorBase::DOFVectorBase(const FiniteElemSpace *f, std::string n) : feSpace(f), name(n), boundaryManager(NULL) { nBasFcts = feSpace->getBasisFcts()->getNumber(); dim = feSpace->getMesh()->getDim(); localIndices.resize(omp_get_overall_max_threads()); localVectors.resize(omp_get_overall_max_threads()); grdPhis.resize(omp_get_overall_max_threads()); grdTmp.resize(omp_get_overall_max_threads()); D2Phis.resize(omp_get_overall_max_threads()); for (int i = 0; i < omp_get_overall_max_threads(); i++) { localIndices[i] = GET_MEMORY(DegreeOfFreedom, this->nBasFcts); localVectors[i] = GET_MEMORY(T, this->nBasFcts); grdPhis[i] = NEW DimVec(dim, DEFAULT_VALUE, 0.0); grdTmp[i] = NEW DimVec(dim, DEFAULT_VALUE, 0.0); D2Phis[i] = NEW DimMat(dim, NO_INIT); } elementVector = NEW ElementVector(this->nBasFcts); } template DOFVectorBase::~DOFVectorBase() { for (int i = 0; i < static_cast(localIndices.size()); i++) { FREE_MEMORY(localIndices[i], DegreeOfFreedom, this->nBasFcts); FREE_MEMORY(localVectors[i], T, this->nBasFcts); DELETE grdPhis[i]; DELETE grdTmp[i]; DELETE D2Phis[i]; } } template DOFVector::DOFVector(const FiniteElemSpace* f, std::string n) : DOFVectorBase(f, n), coarsenOperation(NO_OPERATION) { init(f, n); } template void DOFVector::init(const FiniteElemSpace* f, std::string n) { name = n; feSpace = f; if (feSpace && feSpace->getAdmin()) { (feSpace->getAdmin())->addDOFIndexed(this); } this->boundaryManager = NEW BoundaryManager(f); } template DOFVector::~DOFVector() { if (feSpace && feSpace->getAdmin()) { (feSpace->getAdmin())->removeDOFIndexed(this); } if (this->boundaryManager) { DELETE this->boundaryManager; } vec.clear(); } template DOFVector * DOFVector::traverseVector = NULL; template FastQuadrature *DOFVector::quad_fast = NULL; template double DOFVector::norm = 0.0; template void DOFVectorBase::addElementVector(T factor, const ElementVector &elVec, const BoundaryType *bound, bool add) { FUNCNAME("DOFVector::addElementVector()"); int n_row = elVec.getSize(); for (DegreeOfFreedom i = 0; i < n_row; i++) { BoundaryCondition *condition = bound ? this->getBoundaryManager()->getBoundaryCondition(bound[i]) : NULL; if (!(condition && condition->isDirichlet())) { DegreeOfFreedom irow = elVec.dofIndices[i]; (*this)[irow] = (add ? (*this)[irow] : 0.0); (*this)[irow] += factor * elVec[i]; } } } template double DOFVector::nrm2() const { FUNCNAME("DOFVector::nrm2()"); TEST_EXIT_DBG(feSpace)("feSpace is NULL\n"); TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n"); TEST_EXIT_DBG(static_cast(vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); double nrm = 0.0; Iterator vecIterator(dynamic_cast*>(const_cast*>(this)), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) nrm += (*vecIterator) * (*vecIterator); return(sqrt(nrm)); } template double DOFVector::squareNrm2() const { FUNCNAME("DOFVector::nrm2()"); TEST_EXIT_DBG(feSpace)("feSpace is NULL\n"); TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n"); TEST_EXIT_DBG(static_cast(vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); double nrm = 0.0; Iterator vecIterator(dynamic_cast*>(const_cast*>(this)), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) nrm += (*vecIterator) * (*vecIterator); return nrm; } template T DOFVector::asum() const { FUNCNAME("DOFVector::asum()"); TEST_EXIT_DBG(feSpace)("feSpace is NULL\n"); TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n"); TEST_EXIT_DBG(static_cast(vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); double nrm = 0.0; Iterator vecIterator(dynamic_cast*>(const_cast*>(this)), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) nrm += abs(*vecIterator); return(nrm); } template T DOFVector::sum() const { FUNCNAME("DOFVector::sum()"); TEST_EXIT_DBG(feSpace)("feSpace is NULL\n"); TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n"); TEST_EXIT_DBG(static_cast(vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); double nrm = 0.0; Iterator vecIterator(dynamic_cast*>(const_cast*>(this)), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) nrm += *vecIterator; return(nrm); } template void DOFVector::set(T alpha) { FUNCNAME("DOFVector::set()"); TEST_EXIT_DBG(feSpace)("feSpace is NULL\n"); TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n"); TEST_EXIT_DBG(static_cast(vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); Iterator vecIterator(dynamic_cast*>(this), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) { *vecIterator = alpha ; } } template void DOFVector::copy(const DOFVector& x) { FUNCNAME("DOFVector::copy()"); TEST_EXIT_DBG(feSpace && x.feSpace) ("feSpace is NULL: %8X, %8X\n", feSpace, x.feSpace); TEST_EXIT_DBG(feSpace->getAdmin() && (feSpace->getAdmin() == x.feSpace->getAdmin())) ("no admin or different admins: %8X, %8X\n", feSpace->getAdmin(), x.feSpace->getAdmin()); TEST_EXIT_DBG(static_cast(vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast(x.vec.size()) >= feSpace->getAdmin()->getUsedSize()) ("x.size = %d too small: admin->sizeUsed = %d\n", x.vec.size(), feSpace->getAdmin()->getUsedSize()); Iterator vecIterator(dynamic_cast*>(this), USED_DOFS); Iterator xIterator(dynamic_cast*>(const_cast*>(&x)), USED_DOFS); for (vecIterator.reset(), xIterator.reset(); !vecIterator.end(); ++vecIterator, ++xIterator) { *vecIterator = *xIterator; } } template T DOFVector::min() const { FUNCNAME("DOFVector::min()"); TEST_EXIT_DBG(feSpace && feSpace->getAdmin()) ("pointer is NULL: %8X, %8X\n", this, feSpace->getAdmin()); TEST_EXIT_DBG((static_cast(vec.size())) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); T m; Iterator vecIterator(const_cast*>(dynamic_cast*>(this)), USED_DOFS); for (vecIterator.reset(), m = *vecIterator; !vecIterator.end(); ++vecIterator) { m = std::min(m, *vecIterator); } return m; } template T DOFVector::max() const { FUNCNAME("DOFVector::max()"); TEST_EXIT_DBG(feSpace && feSpace->getAdmin()) ("pointer is NULL: %8X, %8X\n", this, feSpace->getAdmin()); TEST_EXIT_DBG((static_cast(vec.size())) >= feSpace->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), feSpace->getAdmin()->getUsedSize()); T m; Iterator vecIterator(const_cast*>(dynamic_cast*>(this)), USED_DOFS); for (vecIterator.reset(), m = *vecIterator; !vecIterator.end(); ++vecIterator) { m = std::max(m, *vecIterator); } return m; } template T DOFVector::absMax() const { return std::max(abs(max()), abs(min())); } template void gemv(MatrixTranspose transpose, T alpha, const DOFMatrix& a, const DOFVector& x, T beta, DOFVector& y) { FUNCNAME("DOFVector::gemv()"); int j, jcol, ysize; T sum, ax; const DOFMatrix::MatrixRow *row; TEST_EXIT_DBG(a.getRowFESpace() && a.getColFESpace() && x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X, %8X, %8X\n", a.getRowFESpace(), a.getColFESpace(), x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(a.getRowFESpace()->getAdmin() && a.getColFESpace()->getAdmin() && (((transpose == NoTranspose) && (a.getColFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && (a.getRowFESpace()->getAdmin() == y.getFESpace()->getAdmin()))|| ((transpose == Transpose) && (a.getRowFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && (a.getColFESpace()->getAdmin() == y.getFESpace()->getAdmin())))) ("no admin or different admins: %8X, %8X, %8X, %8X\n", a.getRowFESpace()->getAdmin(), a.getColFESpace()->getAdmin(), x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); if (transpose == NoTranspose) { TEST_EXIT_DBG(static_cast(x.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize()) ("x.size = %d too small: admin->sizeUsed = %d\n", x.getSize(), a.getColFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast(y.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("y.size = %d too small: admin->sizeUsed = %d\n", y.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("a.size = %d too small: admin->sizeUsed = %d\n", a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); } else if (transpose == Transpose) { TEST_EXIT_DBG(static_cast(x.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("x.size = %d too small: admin->sizeUsed = %d\n", x.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast(y.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize()) ("y.size = %d too small: admin->sizeUsed = %d\n", y.getSize(), a.getColFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("a.size = %d too small: admin->sizeUsed = %d\n", a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); } ysize = y.getSize(); typename DOFVector::Iterator vecIterator(dynamic_cast*>(&y), FREE_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) { *vecIterator = 0; } if (transpose == NoTranspose) { typename DOFVector::Iterator vecIterator(dynamic_cast*>(&y), USED_DOFS); DOFMatrix::Iterator rowIterator(const_cast(&a), USED_DOFS); for (vecIterator.reset(), rowIterator.reset(); !rowIterator.end(); ++rowIterator, ++vecIterator) { sum = 0; row = &(a[rowIterator.getDOFIndex()]); for (std::vector::iterator colIterator = rowIterator->begin(); colIterator != rowIterator->end(); colIterator++) { jcol = colIterator->col; if (jcol >= 0) { // entry used? sum += (static_cast(colIterator->entry)) * x[jcol]; } else { if (jcol == DOFMatrix::NO_MORE_ENTRIES) break; } } *vecIterator *= beta; *vecIterator += alpha * sum; } } else if (transpose == Transpose) { typename DOFVector::Iterator vecIterator(dynamic_cast*>(&y), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) { *vecIterator *= beta ; } typename DOFVector::Iterator xIterator(dynamic_cast*>(const_cast*>(&x)), USED_DOFS); DOFMatrix::Iterator rowIterator(const_cast(&a), USED_DOFS); for (xIterator.reset(), rowIterator.reset(); !rowIterator.end(); ++rowIterator, ++xIterator) { ax = alpha * (*xIterator); row = &(a[rowIterator.getDOFIndex()]); for (std::vector::iterator colIterator = rowIterator->begin(); colIterator != rowIterator->end(); colIterator++) { jcol = colIterator->col; if (jcol >= 0) // entry used? y[jcol] += ax * (static_cast(colIterator->entry)); else if (jcol == DOFMatrix::NO_MORE_ENTRIES) break; } } } else { ERROR_EXIT("transpose=%d\n", transpose); } } template void DOFVector::print() const { FUNCNAME("DOFVector::print()"); const DOFAdmin *admin = NULL; const char *format; if (feSpace) admin = feSpace->getAdmin(); MSG("Vec `%s':\n", name.c_str()); int j = 0; if (admin) { if (admin->getUsedSize() > 100) format = "%s(%3d,%10.5le)"; else if (admin->getUsedSize() > 10) format = "%s(%2d,%10.5le)"; else format = "%s(%1d,%10.5le)"; Iterator vecIterator(dynamic_cast*>(const_cast*>(this)), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) { if ((j % 3) == 0) { if (j) Msg::print("\n"); MSG(format, "", vecIterator.getDOFIndex(), *vecIterator); } else { Msg::print(format, " ", vecIterator.getDOFIndex(), *vecIterator); } j++; } Msg::print("\n"); } else { MSG("no DOFAdmin, print whole vector.\n"); for (int i = 0; i < static_cast( vec.size()); i++) { if ((j % 3) == 0) { if (j) Msg::print("\n"); MSG("(%d,%10.5e)",i,vec[i]); } else { Msg::print(" (%d,%10.5e)",i,vec[i]); } j++; } Msg::print("\n"); } return; } template int DOFVector::calcMemoryUsage() const { int result = 0; result += sizeof(DOFVector); result += vec.size() * sizeof(T); return result; } template T DOFVectorBase::evalUh(const DimVec& lambda, DegreeOfFreedom* dof_indices) { BasisFunction* phi = const_cast(this->getFESpace()->getBasisFcts()); int numberOfBasFcts = phi->getNumber(); T val = 0.0; for (int i = 0; i < numberOfBasFcts; i++) val += (*this)[dof_indices[i]]*(*phi->getPhi(i))(lambda); return val; } template void DOFVector::interpol(AbstractFunction > *fct) { FUNCNAME("DOFVector::interpol()"); TEST_EXIT_DBG(fct)("no function to interpolate\n"); interFct = fct; if (!this->getFESpace()) { MSG("no dof admin in vec %s, skipping interpolation\n", this->getName().c_str()); return; } if (!(this->getFESpace()->getAdmin())) { MSG("no dof admin in feSpace %s, skipping interpolation\n", this->getFESpace()->getName().c_str()); return; } if (!(this->getFESpace()->getBasisFcts())) { MSG("no basis functions in admin of vec %s, skipping interpolation\n", this->getName().c_str()); return; } if (!(fct)) { MSG("function that should be interpolated only pointer to NULL, "); Msg::print("skipping interpolation\n"); return; } traverseVector = this; TraverseStack stack; ElInfo *elInfo = stack.traverseFirst(this->getFESpace()->getMesh(), -1, Mesh::CALL_LEAF_EL | Mesh::FILL_COORDS); while (elInfo) { interpolFct(elInfo); elInfo = stack.traverseNext(elInfo); } } template int DOFVector::interpolFct(ElInfo* elinfo) { const BasisFunction *basFct = traverseVector->getFESpace()->getBasisFcts(); const DOFAdmin* admin = traverseVector->getFESpace()->getAdmin(); const DegreeOfFreedom *dof = basFct->getLocalIndices(const_cast(elinfo->getElement()), admin, NULL); const T *inter_val = const_cast(basFct)->interpol(elinfo, 0, NULL, traverseVector->interFct, NULL); int number = basFct->getNumber(); for (int i = 0; i < number; i++) (*traverseVector)[dof[i]] = inter_val[i]; return 0; } template double DOFVector::Int(Quadrature* q) const { FUNCNAME("DOFVector::Int()"); Mesh* mesh = feSpace->getMesh(); if (!q) { int deg = 2 * feSpace->getBasisFcts()->getDegree(); q = Quadrature::provideQuadrature(this->dim, deg); } quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), *q, INIT_PHI); norm = 0.0; traverseVector = const_cast*>(this); mesh->traverse(-1, Mesh::CALL_LEAF_EL| Mesh::FILL_COORDS | Mesh::FILL_DET, Int_fct); return norm; } template int DOFVector::Int_fct(ElInfo *elinfo) { double det, normT; const T *uh_vec; int iq; det = elinfo->getDet(); uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL); int numPoints = quad_fast->getNumPoints(); for (normT = iq = 0; iq < numPoints; iq++) { normT += quad_fast->getWeight(iq)*(uh_vec[iq]); } norm += det*normT; return 0; } template double DOFVector::L1Norm(Quadrature* q) const { FUNCNAME("DOFVector::L1Norm()"); if (!q) { int deg = 2 * feSpace->getBasisFcts()->getDegree(); q = Quadrature::provideQuadrature(this->dim, deg); } quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(),*q,INIT_PHI); norm = 0.0; traverseVector = const_cast*>(this); Mesh* mesh = feSpace->getMesh(); mesh->traverse(-1, Mesh::CALL_LEAF_EL| Mesh::FILL_COORDS | Mesh::FILL_DET, L1Norm_fct); return norm; } template int DOFVector::L1Norm_fct(ElInfo *elinfo) { double det, normT; const T *uh_loc, *uh_vec; int iq; det = elinfo->getDet(); uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL); int numPoints = quad_fast->getNumPoints(); for (normT = iq = 0; iq < numPoints; iq++) { normT += quad_fast->getWeight(iq)*abs(uh_vec[iq]); } norm += det*normT; return 0; } template double DOFVector::L2NormSquare(Quadrature* q) const { FUNCNAME("DOFVector::L2NormSquare()"); if (!q) { int deg = 2 * feSpace->getBasisFcts()->getDegree(); q = Quadrature::provideQuadrature(this->dim, deg); } quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), *q, INIT_PHI); norm = 0.0; traverseVector = const_cast*>(this); Mesh* mesh = feSpace->getMesh(); mesh->traverse(-1, Mesh::CALL_LEAF_EL | Mesh::FILL_COORDS | Mesh::FILL_DET, L2Norm_fct); return norm; } template int DOFVector::L2Norm_fct(ElInfo *elinfo) { double det, normT; const T *uh_vec; int iq; det = elinfo->getDet(); uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL); int numPoints = quad_fast->getNumPoints(); for (normT = iq = 0; iq < numPoints; iq++) { normT += quad_fast->getWeight(iq)*sqr(uh_vec[iq]); } norm += det*normT; return 0; } template int DOFVector::H1Norm_fct(ElInfo *elinfo) { double norm2, normT; const WorldVector *grduh_vec; int iq, j; double det = elinfo->getDet(); grduh_vec = traverseVector->getGrdAtQPs(elinfo, NULL, quad_fast, NULL); int dimOfWorld = Global::getGeo(WORLD); int numPoints = quad_fast->getNumPoints(); for (normT = iq = 0; iq < numPoints; iq++) { for (norm2 = j = 0; j < dimOfWorld; j++) norm2 += sqr(grduh_vec[iq][j]); normT += quad_fast->getWeight(iq)*norm2; } norm += det*normT; return 0; } template double DOFVector::H1NormSquare(Quadrature *q) const { FUNCNAME("DOFVector::H1NormSquare()"); if (!q) { int deg = 2 * feSpace->getBasisFcts()->getDegree() - 2; q = Quadrature::provideQuadrature(this->dim, deg); } quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), *q, INIT_GRD_PHI); norm = 0.0; traverseVector = const_cast*>(this); Mesh *mesh = feSpace->getMesh(); mesh->traverse(-1, Mesh::CALL_LEAF_EL | Mesh::FILL_COORDS | Mesh::FILL_DET | Mesh::FILL_GRD_LAMBDA, H1Norm_fct); return norm; } template void DOFVector::compressDOFIndexed(int first, int last, std::vector &newDOF) { for (int i = first; i <= last; i++) { if (newDOF[i] >= 0) { vec[newDOF[i]] = vec[i]; } } } template Flag DOFVectorBase::getAssembleFlag() { Flag fillFlag(0); std::vector::iterator op; for (op = this->operators.begin(); op != this->operators.end(); ++op) { fillFlag |= (*op)->getFillFlag(); } return fillFlag; } template DOFVector& DOFVector::operator=(const DOFVector& rhs ) { feSpace = rhs.feSpace; DOFVectorBase::feSpace = rhs.feSpace; this->nBasFcts = rhs.nBasFcts; vec = rhs.vec; this->elementVector = new ElementVector(this->nBasFcts); interFct = rhs.interFct; coarsenOperation = rhs.coarsenOperation; this->operators = rhs.operators; this->operatorFactor = rhs.operatorFactor; if (rhs.boundaryManager) { if (this->boundaryManager) delete this->boundaryManager; this->boundaryManager = new BoundaryManager(*rhs.boundaryManager); // boundaryManager->setDOFVector(this); } else { this->boundaryManager = NULL; } return *this; } template const DOFVector& operator*=(DOFVector& x, T scal) { FUNCNAME("DOFVector::operator*=(DOFVector& x, T scal)"); TEST_EXIT_DBG(x.getFESpace() && x.getFESpace()->getAdmin()) ("pointer is NULL: %8X, %8X\n", x.getFESpace(), x.getFESpace()->getAdmin()); typename DOFVector::Iterator vecIterator(dynamic_cast*>(&x), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) { (*vecIterator) *= scal; } return x; } template const DOFVector& operator+=(DOFVector& x, const DOFVector& y) { FUNCNAME("DOFVector::operator+=(DOFVector& x, const DOFVector& y)"); TEST_EXIT_DBG(x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(x.getFESpace()->getAdmin() && (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin())) ("no admin or different admins: %8X, %8X\n", x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n"); typename DOFVector::Iterator xIterator(dynamic_cast*>(&x), USED_DOFS); typename DOFVector::Iterator yIterator(dynamic_cast*>(const_cast*>(&y)), USED_DOFS); for (xIterator.reset(), yIterator.reset(); !xIterator.end(); ++xIterator, ++yIterator) { *xIterator += *yIterator; } return x; } template const DOFVector& operator-=(DOFVector& x, const DOFVector& y) { FUNCNAME("DOFVector::operator-=(DOFVector& x, const DOFVector& y)"); TEST_EXIT_DBG(x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(x.getFESpace()->getAdmin() && (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin())) ("no admin or different admins: %8X, %8X\n", x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n"); typename DOFVector::Iterator xIterator(dynamic_cast*>(&x), USED_DOFS); typename DOFVector::Iterator yIterator(dynamic_cast*>(const_cast*>(&y)), USED_DOFS); for (xIterator.reset(), yIterator.reset(); !xIterator.end(); ++xIterator, ++yIterator) { *xIterator -= *yIterator; } return x; } template const DOFVector& operator*=(DOFVector& x, const DOFVector& y) { FUNCNAME("DOFVector::operator*=(DOFVector& x, const DOFVector& y)"); TEST_EXIT_DBG(x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(x.getFESpace()->getAdmin() && (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin())) ("no admin or different admins: %8X, %8X\n", x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n"); typename DOFVector::Iterator xIterator(dynamic_cast*>(&x), USED_DOFS); typename DOFVector::Iterator yIterator(dynamic_cast*>(const_cast*>(&y)), USED_DOFS); for (xIterator.reset(), yIterator.reset(); !xIterator.end(); ++xIterator, ++yIterator) { *xIterator *= *yIterator; } return x; } template T operator*(DOFVector& x, DOFVector& y) { FUNCNAME("DOFVector::operator*(DOFVector& x, DOFVector& y)"); TEST_EXIT_DBG(x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(x.getFESpace()->getAdmin() && (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin())) ("no admin or different admins: %8X, %8X\n", x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n"); T dot = 0; typename DOFVector::Iterator xIterator(dynamic_cast*>(&x), USED_DOFS); typename DOFVector::Iterator yIterator(dynamic_cast*>(&y), USED_DOFS); for (xIterator.reset(), yIterator.reset(); !xIterator.end(); ++xIterator, ++yIterator) { dot += (*xIterator) * (*yIterator); } return(dot); } template void mv(MatrixTranspose transpose, const DOFMatrix &a, const DOFVector&x, DOFVector &result, bool add) { FUNCNAME("DOFVector::mv()"); TEST_EXIT_DBG(a.getRowFESpace() && a.getColFESpace() && x.getFESpace() && result.getFESpace()) ("getFESpace() is NULL: %8X, %8X, %8X, %8X\n", a.getRowFESpace(), a.getColFESpace(), x.getFESpace(), result.getFESpace()); TEST_EXIT_DBG((a.getRowFESpace()->getAdmin() && a.getColFESpace()->getAdmin()) && (((transpose == NoTranspose) && (a.getColFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && (a.getRowFESpace()->getAdmin() == result.getFESpace()->getAdmin())) || ((transpose == Transpose) && (a.getRowFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && (a.getColFESpace()->getAdmin() == result.getFESpace()->getAdmin())))) ("no admin or different admins: %8X, %8X, %8X, %8X\n", a.getRowFESpace()->getAdmin(), a.getColFESpace()->getAdmin(), x.getFESpace()->getAdmin(), result.getFESpace()->getAdmin()); if (transpose == NoTranspose) { TEST_EXIT_DBG(static_cast(x.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize()) ("x.size = %d too small: admin->sizeUsed = %d\n", x.getSize(), a.getColFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast( result.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", result.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("a.size = %d too small: admin->sizeUsed = %d\n", a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); typename DOFVector::Iterator vecIterator(dynamic_cast*>(&result), USED_DOFS); DOFMatrix::Iterator rowIterator(const_cast(&a), USED_DOFS); for (vecIterator.reset(), rowIterator.reset(); !rowIterator.end(); ++rowIterator, ++vecIterator) { double sum = 0; if (!add) *vecIterator = 0.0; for (std::vector::iterator colIterator = rowIterator->begin(); colIterator != rowIterator->end(); colIterator++) { int jcol = colIterator->col; if (jcol >= 0) { // entry used? sum += (static_cast(colIterator->entry)) * x[jcol]; } else { if (jcol == DOFMatrix::NO_MORE_ENTRIES) break; } } *vecIterator += sum; } } else if (transpose == Transpose) { TEST_EXIT_DBG(static_cast(x.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("x.size = %d too small: admin->sizeUsed = %d\n", x.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast( result.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize()) ("size = %d too small: admin->sizeUsed = %d\n", result.getSize(), a.getColFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize()) ("a.size = %d too small: admin->sizeUsed = %d\n", a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize()); if (!add) { typename DOFVector::Iterator vecIterator(dynamic_cast*>(&result), USED_DOFS); for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) { *vecIterator = 0; } } typename DOFVector::Iterator xIterator(const_cast*>(dynamic_cast*>(&x)), USED_DOFS); DOFMatrix::Iterator rowIterator(const_cast(&a), USED_DOFS); for (xIterator.reset(), rowIterator.reset(); !rowIterator.end(); ++rowIterator, ++xIterator) { T ax = (*xIterator); for (std::vector::iterator colIterator = rowIterator->begin(); colIterator != rowIterator->end(); colIterator++) { int jcol = colIterator->col; if (jcol >= 0) // entry used? result[jcol] += (static_cast(colIterator->entry)) * ax; else if (jcol == DOFMatrix::NO_MORE_ENTRIES) break; } } } else { ERROR_EXIT("transpose = %d\n", transpose); } } template void axpy(double alpha, const DOFVector& x, DOFVector& y) { FUNCNAME("DOFVector::axpy()"); TEST_EXIT_DBG(x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(x.getFESpace()->getAdmin() && (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin())) ("no admin or different admins: %8X, %8X\n", x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); TEST_EXIT_DBG(static_cast(x.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize()) ("size = %d too small: admin->size = %d\n", x.getSize(), x.getFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast(y.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize()) ("y.size = %d too small: admin->size = %d\n", y.getSize(), x.getFESpace()->getAdmin()->getUsedSize()); typename DOFVector::Iterator xIterator(dynamic_cast*>(const_cast*>(&x)), USED_DOFS); typename DOFVector::Iterator yIterator(dynamic_cast*>(&y), USED_DOFS); for (xIterator.reset(), yIterator.reset(); !xIterator.end(); ++xIterator, ++yIterator) { *yIterator += alpha * (*xIterator); } } template const DOFVector& operator*(const DOFVector& v, double d) { static DOFVector result; return mult(d, v, result); } template const DOFVector& operator*(double d, const DOFVector& v) { static DOFVector result; return mult(d, v, result); } template const DOFVector& operator+(const DOFVector &v1 , const DOFVector &v2) { static DOFVector result; return add(v1, v2, result); } template void xpay(double alpha, const DOFVector& x, DOFVector& y) { FUNCNAME("DOFVector::xpay()"); TEST_EXIT_DBG(x.getFESpace() && y.getFESpace()) ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace()); TEST_EXIT_DBG(x.getFESpace()->getAdmin() && (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin())) ("no admin or different admins: %8X, %8X\n", x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin()); TEST_EXIT_DBG(static_cast(x.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize()) ("size = %d too small: admin->size = %d\n", x.getSize(), x.getFESpace()->getAdmin()->getUsedSize()); TEST_EXIT_DBG(static_cast(y.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize()) ("y.size = %d too small: admin->size = %d\n", y.getSize(), x.getFESpace()->getAdmin()->getUsedSize()); typename DOFVector::Iterator xIterator(dynamic_cast*>(const_cast*>(&x)), USED_DOFS); typename DOFVector::Iterator yIterator(dynamic_cast*>(&y), USED_DOFS); for (xIterator.reset(), yIterator.reset(); !xIterator.end(); ++xIterator, ++yIterator) { *yIterator = alpha *(*yIterator)+ (*xIterator); } } template inline const DOFVector& mult(double scal, const DOFVector& v, DOFVector& result) { typename DOFVector::Iterator vIterator(dynamic_cast*>(const_cast*>(&v)), USED_DOFS); typename DOFVector::Iterator rIterator(dynamic_cast*>(&result), USED_DOFS); for (vIterator.reset(), rIterator.reset(); !vIterator.end(); ++vIterator, ++rIterator) { *rIterator = scal * (*vIterator); } return result; } template inline const DOFVector& add(const DOFVector& v, double scal, DOFVector& result) { typename DOFVector::Iterator vIterator(dynamic_cast*>(const_cast*>(&v)), USED_DOFS); typename DOFVector::Iterator rIterator(dynamic_cast*>(&result), USED_DOFS); for (vIterator.reset(), rIterator.reset(); !vIterator.end(); ++vIterator, ++rIterator) { *rIterator = (*vIterator) + scal; } return result; } template inline const DOFVector& add(const DOFVector& v1, const DOFVector& v2, DOFVector& result) { typename DOFVector::Iterator v1Iterator(dynamic_cast*>(const_cast*>(&v1)), USED_DOFS); typename DOFVector::Iterator v2Iterator(dynamic_cast*>(const_cast*>(&v2)), USED_DOFS); typename DOFVector::Iterator rIterator(dynamic_cast*>(&result), USED_DOFS); for (v1Iterator.reset(), v2Iterator.reset(), rIterator.reset(); !v1Iterator.end(); ++v1Iterator, ++v2Iterator, ++rIterator) { *rIterator = (*v1Iterator) + (*v2Iterator); } return result; } template inline const DOFVector& mod(const DOFVector& v, DOFVector& result) { typename DOFVector::Iterator vIterator(dynamic_cast*>(const_cast*>(&v)), USED_DOFS); typename DOFVector::Iterator rIterator(dynamic_cast*>(&result), USED_DOFS); for (vIterator.reset(), rIterator.reset(); !vIterator.end(); ++vIterator, ++rIterator) { *rIterator = fmod((*vIterator), 1.0); } return result; } template inline const DOFVector& Tanh(const DOFVector& v, DOFVector& result) { typename DOFVector::Iterator vIterator(dynamic_cast*>(const_cast*>(&v)), USED_DOFS); typename DOFVector::Iterator rIterator(dynamic_cast*>(&result), USED_DOFS); double eps; GET_PARAMETER(1, "vecphase->epsilon", "%f", &eps); for (vIterator.reset(), rIterator.reset(); !vIterator.end(); ++vIterator, ++rIterator) { *rIterator = 0.5 * (1.0 - tanh(3.0 * (*vIterator) / eps)); } return result; } template const T *DOFVectorBase::getVecAtQPs(const ElInfo *elInfo, const Quadrature *quad, const FastQuadrature *quadFast, T *vecAtQPs) const { FUNCNAME("DOFVector::getVecAtQPs()"); TEST_EXIT_DBG(quad || quadFast)("neither quad nor quadFast defined\n"); if (quad && quadFast) { TEST_EXIT_DBG(quad == quadFast->getQuadrature()) ("quad != quadFast->quadrature\n"); } TEST_EXIT_DBG(!quadFast || quadFast->getBasisFunctions() == feSpace->getBasisFcts()) ("invalid basis functions"); const BasisFunction *basFcts = feSpace->getBasisFcts(); int nPoints = quadFast ? quadFast->getQuadrature()->getNumPoints() : quad->getNumPoints(); static T *localvec = NULL; T *result; if (vecAtQPs) { result = vecAtQPs; } else { if (localvec) delete [] localvec; localvec = new T[nPoints]; for (int i = 0; i < nPoints; i++) { localvec[i] = 0.0; } result = localvec; } T *localVec = localVectors[omp_get_thread_num()]; getLocalVector(elInfo->getElement(), localVec); for (int i = 0; i < nPoints; i++) { result[i] = 0.0; for (int j = 0; j < nBasFcts; j++) { result[i] += localVec[j] * (quadFast ? (quadFast->getPhi(i, j)) : ((*(basFcts->getPhi(j)))(quad->getLambda(i)))); } } return const_cast(result); } template const T *DOFVectorBase::getLocalVector(const Element *el, T *d) const { static T* localVec = NULL; T *result; if (d) { result = d; } else { if (localVec) delete [] localVec; localVec = new T[nBasFcts]; result = localVec; } DegreeOfFreedom *myLocalIndices = localIndices[omp_get_thread_num()]; feSpace->getBasisFcts()->getLocalIndices(el, feSpace->getAdmin(), myLocalIndices); for (int i = 0; i < nBasFcts; i++) { result[i] = (*this)[myLocalIndices[i]]; } return result; } // integral u^2(1-u)^2 template double DOFVector::DoubleWell(Quadrature* q) const { FUNCNAME("DOFVector::DoubleWell()"); if (!q) { int deg = 2 * feSpace->getBasisFcts()->getDegree(); q = Quadrature::provideQuadrature(this->dim, deg); } quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), *q, INIT_PHI); norm = 0.0; traverseVector = const_cast*>(this); Mesh* mesh = feSpace->getMesh(); mesh->traverse(-1, Mesh::CALL_LEAF_EL | Mesh::FILL_COORDS | Mesh::FILL_DET, DoubleWell_fct); return norm; } template int DOFVector::DoubleWell_fct(ElInfo *elinfo) { const T *uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL); int nPoints = quad_fast->getNumPoints(); double normT = 0.0; for (int iq = 0; iq < nPoints; iq++) { normT += quad_fast->getWeight(iq) * sqr(uh_vec[iq]) * sqr(1.0 - uh_vec[iq]); } norm += elinfo->getDet() * normT; return 0; } }