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; }