bool ASMu2D::globalL2projection (Matrix& sField, const IntegrandBase& integrand, bool continuous) const { if (!lrspline) return true; // silently ignore empty patches PROFILE2("ASMu2D::globalL2"); const int p1 = lrspline->order(0); const int p2 = lrspline->order(1); // Get Gaussian quadrature points const int ng1 = continuous ? nGauss : p1 - 1; const int ng2 = continuous ? nGauss : p2 - 1; const double* xg = GaussQuadrature::getCoord(ng1); const double* yg = GaussQuadrature::getCoord(ng2); const double* wg = continuous ? GaussQuadrature::getWeight(nGauss) : nullptr; if (!xg || !yg) return false; if (continuous && !wg) return false; // Set up the projection matrices const size_t nnod = this->getNoNodes(); const size_t ncomp = integrand.getNoFields(); SparseMatrix A(SparseMatrix::SUPERLU); StdVector B(nnod*ncomp); A.redim(nnod,nnod); double dA = 0.0; Vector phi; Matrix dNdu, Xnod, Jac; Go::BasisDerivsSf spl1; Go::BasisPtsSf spl0; // === Assembly loop over all elements in the patch ========================== std::vector<LR::Element*>::iterator el = lrspline->elementBegin(); for (int iel = 1; el != lrspline->elementEnd(); el++, iel++) { if (continuous) { // Set up control point (nodal) coordinates for current element if (!this->getElementCoordinates(Xnod,iel)) return false; else if ((dA = 0.25*this->getParametricArea(iel)) < 0.0) return false; // topology error (probably logic error) } // Compute parameter values of the Gauss points over this element std::array<RealArray,2> gpar, unstrGpar; this->getGaussPointParameters(gpar[0],0,ng1,iel,xg); this->getGaussPointParameters(gpar[1],1,ng2,iel,yg); // convert to unstructred mesh representation expandTensorGrid(gpar.data(),unstrGpar.data()); // Evaluate the secondary solution at all integration points if (!this->evalSolution(sField,integrand,unstrGpar.data())) return false; // set up basis function size (for extractBasis subroutine) phi.resize((**el).nBasisFunctions()); // --- Integration loop over all Gauss points in each direction ---------- int ip = 0; for (int j = 0; j < ng2; j++) for (int i = 0; i < ng1; i++, ip++) { if (continuous) { lrspline->computeBasis(gpar[0][i],gpar[1][j],spl1,iel-1); SplineUtils::extractBasis(spl1,phi,dNdu); } else { lrspline->computeBasis(gpar[0][i],gpar[1][j],spl0,iel-1); phi = spl0.basisValues; } // Compute the Jacobian inverse and derivatives double dJw = 1.0; if (continuous) { dJw = dA*wg[i]*wg[j]*utl::Jacobian(Jac,dNdu,Xnod,dNdu,false); if (dJw == 0.0) continue; // skip singular points } // Integrate the linear system A*x=B for (size_t ii = 0; ii < phi.size(); ii++) { int inod = MNPC[iel-1][ii]+1; for (size_t jj = 0; jj < phi.size(); jj++) { int jnod = MNPC[iel-1][jj]+1; A(inod,jnod) += phi[ii]*phi[jj]*dJw; } for (size_t r = 1; r <= ncomp; r++) B(inod+(r-1)*nnod) += phi[ii]*sField(r,ip+1)*dJw; } } } #if SP_DEBUG > 2 std::cout <<"---- Matrix A -----"<< A <<"-------------------"<< std::endl; std::cout <<"---- Vector B -----"<< B <<"-------------------"<< std::endl; #endif // Solve the patch-global equation system if (!A.solve(B)) return false; // Store the control-point values of the projected field sField.resize(ncomp,nnod); for (size_t i = 1; i <= nnod; i++) for (size_t j = 1; j <= ncomp; j++) sField(j,i) = B(i+(j-1)*nnod); return true; }
bool ASMu3Dmx::assembleL2matrices (SparseMatrix& A, StdVector& B, const IntegrandBase& integrand, bool continuous) const { const int p1 = projBasis->order(0); const int p2 = projBasis->order(1); const int p3 = projBasis->order(2); // Get Gaussian quadrature points const int ng1 = continuous ? nGauss : p1 - 1; const int ng2 = continuous ? nGauss : p2 - 1; const int ng3 = continuous ? nGauss : p3 - 1; const double* xg = GaussQuadrature::getCoord(ng1); const double* yg = GaussQuadrature::getCoord(ng2); const double* zg = GaussQuadrature::getCoord(ng3); const double* wg = continuous ? GaussQuadrature::getWeight(nGauss) : nullptr; if (!xg || !yg || !zg) return false; if (continuous && !wg) return false; size_t nnod = this->getNoProjectionNodes(); double dV = 0.0; Vectors phi(2); Matrices dNdu(2); Matrix sField, Xnod, Jac; std::vector<Go::BasisDerivs> spl1(2); std::vector<Go::BasisPts> spl0(2); // === Assembly loop over all elements in the patch ========================== LR::LRSplineVolume* geoVol; if (m_basis[geoBasis-1]->nBasisFunctions() == projBasis->nBasisFunctions()) geoVol = m_basis[geoBasis-1].get(); else geoVol = projBasis.get(); for (const LR::Element* el1 : geoVol->getAllElements()) { double uh = (el1->umin()+el1->umax())/2.0; double vh = (el1->vmin()+el1->vmax())/2.0; double wh = (el1->wmin()+el1->wmax())/2.0; std::vector<size_t> els; els.push_back(projBasis->getElementContaining(uh, vh, wh) + 1); els.push_back(m_basis[geoBasis-1]->getElementContaining(uh, vh, wh) + 1); if (continuous) { // Set up control point (nodal) coordinates for current element if (!this->getElementCoordinates(Xnod,els[1])) return false; else if ((dV = 0.25*this->getParametricVolume(els[1])) < 0.0) return false; // topology error (probably logic error) } // Compute parameter values of the Gauss points over this element RealArray gpar[3], unstrGpar[3]; this->getGaussPointParameters(gpar[0],0,ng1,els[1],xg); this->getGaussPointParameters(gpar[1],1,ng2,els[1],yg); this->getGaussPointParameters(gpar[2],2,ng3,els[1],zg); // convert to unstructred mesh representation expandTensorGrid(gpar, unstrGpar); // Evaluate the secondary solution at all integration points if (!this->evalSolution(sField,integrand,unstrGpar)) return false; // set up basis function size (for extractBasis subroutine) const LR::Element* elm = projBasis->getElement(els[0]-1); phi[0].resize(elm->nBasisFunctions()); phi[1].resize(el1->nBasisFunctions()); IntVec lmnpc; if (projBasis != m_basis[0]) { lmnpc.reserve(phi[0].size()); for (const LR::Basisfunction* f : elm->support()) lmnpc.push_back(f->getId()); } const IntVec& mnpc = projBasis == m_basis[0] ? MNPC[els[1]-1] : lmnpc; // --- Integration loop over all Gauss points in each direction ---------- Matrix eA(phi[0].size(), phi[0].size()); Vectors eB(sField.rows(), Vector(phi[0].size())); int ip = 0; for (int k = 0; k < ng3; k++) for (int j = 0; j < ng2; j++) for (int i = 0; i < ng1; i++, ip++) { if (continuous) { projBasis->computeBasis(gpar[0][i], gpar[1][j], gpar[2][k], spl1[0], els[0]-1); SplineUtils::extractBasis(spl1[0],phi[0],dNdu[0]); m_basis[geoBasis-1]->computeBasis(gpar[0][i], gpar[1][j], gpar[2][k], spl1[1], els[1]-1); SplineUtils::extractBasis(spl1[1], phi[1], dNdu[1]); } else { projBasis->computeBasis(gpar[0][i], gpar[1][j], gpar[2][k], spl0[0], els[0]-1); phi[0] = spl0[0].basisValues; } // Compute the Jacobian inverse and derivatives double dJw = 1.0; if (continuous) { dJw = dV*wg[i]*wg[j]*wg[k]*utl::Jacobian(Jac,dNdu[1],Xnod,dNdu[1],false); if (dJw == 0.0) continue; // skip singular points } // Integrate the mass matrix eA.outer_product(phi[0], phi[0], true, dJw); // Integrate the rhs vector B for (size_t r = 1; r <= sField.rows(); r++) eB[r-1].add(phi[0],sField(r,ip+1)*dJw); } for (size_t i = 0; i < eA.rows(); ++i) { for (size_t j = 0; j < eA.cols(); ++j) A(mnpc[i]+1, mnpc[j]+1) += eA(i+1,j+1); int jp = mnpc[i]+1; for (size_t r = 0; r < sField.rows(); r++, jp += nnod) B(jp) += eB[r](1+i); } } return true; }
LR::LRSplineSurface* ASMu2D::scRecovery (const IntegrandBase& integrand) const { PROFILE2("ASMu2D::scRecovery"); const int m = integrand.derivativeOrder(); const int p1 = lrspline->order(0); const int p2 = lrspline->order(1); // Get Gaussian quadrature point coordinates const int ng1 = p1 - m; const int ng2 = p2 - m; const double* xg = GaussQuadrature::getCoord(ng1); const double* yg = GaussQuadrature::getCoord(ng2); if (!xg || !yg) return nullptr; // Compute parameter values of the Greville points std::array<RealArray,2> gpar; if (!this->getGrevilleParameters(gpar[0],0)) return nullptr; if (!this->getGrevilleParameters(gpar[1],1)) return nullptr; const int n1 = p1 - m + 1; // Patch size in first parameter direction const int n2 = p2 - m + 1; // Patch size in second parameter direction const size_t nCmp = integrand.getNoFields(); // Number of result components const size_t nPol = n1*n2; // Number of terms in polynomial expansion Matrix sValues(nCmp,gpar[0].size()); Vector P(nPol); Go::Point X, G; // Loop over all Greville points (one for each basis function) size_t k, l, ip = 0; std::vector<LR::Element*>::const_iterator elStart, elEnd, el; std::vector<LR::Element*> supportElements; for (LR::Basisfunction *b : lrspline->getAllBasisfunctions()) { #if SP_DEBUG > 2 std::cout <<"Basis: "<< *b <<"\n ng1 ="<< ng1 <<"\n ng2 ="<< ng2 <<"\n nPol="<< nPol << std::endl; #endif // Special case for basis functions with too many zero knot spans by using // the extended support // if(nel*ng1*ng2 < nPol) if(true) { // KMO: Here I'm not sure how this will change when m > 1. // In that case I think we would need smaller patches (as in the tensor // splines case). But how to do that??? supportElements = b->getExtendedSupport(); elStart = supportElements.begin(); elEnd = supportElements.end(); #if SP_DEBUG > 2 std::cout <<"Extended basis:"; for (el = elStart; el != elEnd; el++) std::cout <<"\n " << **el; std::cout << std::endl; #endif } else { elStart = b->supportedElementBegin(); elEnd = b->supportedElementEnd(); } // Physical coordinates of current Greville point lrspline->point(G,gpar[0][ip],gpar[1][ip]); // Set up the local projection matrices DenseMatrix A(nPol,nPol); Matrix B(nPol,nCmp); // Loop over all non-zero knot-spans in the support of // the basis function associated with current Greville point for (el = elStart; el != elEnd; el++) { int iel = (**el).getId()+1; // evaluate all gauss points for this element std::array<RealArray,2> gaussPt, unstrGauss; this->getGaussPointParameters(gaussPt[0],0,ng1,iel,xg); this->getGaussPointParameters(gaussPt[1],1,ng2,iel,yg); #if SP_DEBUG > 2 std::cout << "Element " << **el << std::endl; #endif // convert to unstructred mesh representation expandTensorGrid(gaussPt.data(),unstrGauss.data()); // Evaluate the secondary solution at all Gauss points Matrix sField; if (!this->evalSolution(sField,integrand,unstrGauss.data())) return nullptr; // Loop over the Gauss points in current knot-span int i, j, ig = 1; for (j = 0; j < ng2; j++) for (i = 0; i < ng1; i++, ig++) { // Evaluate the polynomial expansion at current Gauss point lrspline->point(X,gaussPt[0][i],gaussPt[1][j]); evalMonomials(n1,n2,X[0]-G[0],X[1]-G[1],P); #if SP_DEBUG > 2 std::cout <<"Greville point: "<< G <<"\nGauss point: "<< gaussPt[0][i] <<", "<< gaussPt[0][j] <<"\nMapped gauss point: "<< X <<"\nP-matrix:"<< P <<"--------------------\n"<< std::endl; #endif for (k = 1; k <= nPol; k++) { // Accumulate the projection matrix, A += P^t * P for (l = 1; l <= nPol; l++) A(k,l) += P(k)*P(l); // Accumulate the right-hand-side matrix, B += P^t * sigma for (l = 1; l <= nCmp; l++) B(k,l) += P(k)*sField(l,ig); } } } #if SP_DEBUG > 2 std::cout <<"---- Matrix A -----"<< A <<"-------------------"<< std::endl; std::cout <<"---- Vector B -----"<< B <<"-------------------"<< std::endl; #endif // Solve the local equation system if (!A.solve(B)) return nullptr; // Evaluate the projected field at current Greville point (first row of B) ip++; for (l = 1; l <= nCmp; l++) sValues(l,ip) = B(1,l); #if SP_DEBUG > 1 std::cout <<"Greville point "<< ip <<" (u,v) = " << gpar[0][ip-1] <<" "<< gpar[1][ip-1] <<" :"; for (k = 1; k <= nCmp; k++) std::cout <<" "<< sValues(k,ip); std::cout << std::endl; #endif } // Project the Greville point results onto the spline basis // to find the control point values return this->regularInterpolation(gpar[0],gpar[1],sValues); }