void fullMatrix<double>::gemm(const fullMatrix<double> &a, const fullMatrix<double> &b, double alpha, double beta, bool transposeA, bool transposeB) { int M = size1(), N = size2(), K = transposeA ? a.size1() : a.size2(); int LDA = a.size1(), LDB = b.size1(), LDC = size1(); F77NAME(dgemm)(transposeA ? "T" : "N", transposeB ? "T" :"N", &M, &N, &K, &alpha, a._data, &LDA, b._data, &LDB, &beta, _data, &LDC); }
void fullMatrix<double>::multOnBlock(const fullMatrix<double> &b, const int ncol, const int fcol, const int alpha_, const int beta_, fullVector<double> &c) const { int M = 1, N = ncol, K = b.size1() ; int LDA = _r, LDB = b.size1(), LDC = 1; double alpha = alpha_, beta = beta_; F77NAME(dgemm)("N", "N", &M, &N, &K, &alpha, _data, &LDA, &(b._data[fcol*K]), &LDB, &beta, &(c._data[fcol]), &LDC); }
void fullMatrix<double>::mult(const fullMatrix<double> &b, fullMatrix<double> &c) const { int M = c.size1(), N = c.size2(), K = _c; int LDA = _r, LDB = b.size1(), LDC = c.size1(); double alpha = 1., beta = 0.; F77NAME(dgemm)("N", "N", &M, &N, &K, &alpha, _data, &LDA, b._data, &LDB, &beta, c._data, &LDC); }
void PViewFactory::setEntry (int id, const fullMatrix<double> &val) { std::vector<double> &vv = _values[id]; vv.resize(val.size1()*val.size2()); int k=0; for (int i=0;i<val.size1(); i++) { for (int j=0;j<val.size2(); j++) { vv[k++] = val(i,j); } } }
bool fullMatrix<double>::invert(fullMatrix<double> &result) const { int M = size1(), N = size2(), lda = size1(), info; int *ipiv = new int[std::min(M, N)]; if (result.size2() != M || result.size1() != N) { if (result._own_data || !result._data) result.resize(M,N,false); else Msg::Fatal("FullMatrix: Bad dimension, I cannot write in proxy"); } result.setAll(*this); F77NAME(dgetrf)(&M, &N, result._data, &lda, ipiv, &info); if(info == 0){ int lwork = M * 4; double *work = new double[lwork]; F77NAME(dgetri)(&M, result._data, &lda, ipiv, work, &lwork, &info); delete [] work; } delete [] ipiv; if(info == 0) return true; else if(info > 0) Msg::Error("U(%d,%d)=0 in matrix inversion", info, info); else Msg::Error("Wrong %d-th argument in matrix inversion", -info); return false; }
void polynomialBasis::df(const fullMatrix<double> &coord, fullMatrix<double> &dfm) const { double dfv[1256][3]; dfm.resize(coord.size1() * 3, coefficients.size1(), false); int dimCoord = coord.size2(); for(int iPoint = 0; iPoint < coord.size1(); iPoint++) { df(coord(iPoint, 0), dimCoord > 1 ? coord(iPoint, 1) : 0., dimCoord > 2 ? coord(iPoint, 2) : 0., dfv); for(int i = 0; i < coefficients.size1(); i++) { dfm(iPoint * 3 + 0, i) = dfv[i][0]; dfm(iPoint * 3 + 1, i) = dfv[i][1]; dfm(iPoint * 3 + 2, i) = dfv[i][2]; } } }
void polynomialBasis::f(const fullMatrix<double> &coord, fullMatrix<double> &sf) const { double p[1256]; sf.resize(coord.size1(), coefficients.size1()); for(int iPoint = 0; iPoint < coord.size1(); iPoint++) { evaluateMonomials(coord(iPoint, 0), coord.size2() > 1 ? coord(iPoint, 1) : 0, coord.size2() > 2 ? coord(iPoint, 2) : 0, p); for(int i = 0; i < coefficients.size1(); i++) { sf(iPoint, i) = 0.; for(int j = 0; j < coefficients.size2(); j++) sf(iPoint, i) += coefficients(i, j) * p[j]; } } }
void linearSystemPETSc<fullMatrix<double> >::addToMatrix(int row, int col, const fullMatrix<double> &val) { if (!_entriesPreAllocated) preAllocateEntries(); #ifdef PETSC_USE_COMPLEX fullMatrix<std::complex<double> > modval(val.size1(), val.size2()); for (int ii = 0; ii < val.size1(); ii++) { for (int jj = 0; jj < val.size1(); jj++) { modval(ii, jj) = val (jj, ii); modval(jj, ii) = val (ii, jj); } } #else fullMatrix<double> &modval = *const_cast<fullMatrix<double> *>(&val); for (int ii = 0; ii < val.size1(); ii++) { for (int jj = 0; jj < ii; jj++) { PetscScalar buff = modval(ii, jj); modval(ii, jj) = modval (jj, ii); modval(jj, ii) = buff; } } #endif PetscInt i = row, j = col; MatSetValuesBlocked(_a, 1, &i, 1, &j, &modval(0,0), ADD_VALUES); //transpose back so that the original matrix is not modified #ifndef PETSC_USE_COMPLEX for (int ii = 0; ii < val.size1(); ii++) for (int jj = 0; jj < ii; jj++) { PetscScalar buff = modval(ii,jj); modval(ii, jj) = modval (jj,ii); modval(jj, ii) = buff; } #endif _valuesNotAssembled = true; }
bool fullMatrix<double>::svd(fullMatrix<double> &V, fullVector<double> &S) { fullMatrix<double> VT(V.size2(), V.size1()); int M = size1(), N = size2(), LDA = size1(), LDVT = VT.size1(), info; int lwork = std::max(3 * std::min(M, N) + std::max(M, N), 5 * std::min(M, N)); fullVector<double> WORK(lwork); F77NAME(dgesvd)("O", "A", &M, &N, _data, &LDA, S._data, _data, &LDA, VT._data, &LDVT, WORK._data, &lwork, &info); V = VT.transpose(); if(info == 0) return true; if(info > 0) Msg::Error("SVD did not converge"); else Msg::Error("Wrong %d-th argument in SVD decomposition", -info); return false; }
void ana(double (*f)(fullVector<double>& xyz), const fullMatrix<double>& point, fullMatrix<double>& eval){ // Alloc eval for Scalar Values // const size_t nPoint = point.size1(); eval.resize(nPoint, 1); // Loop on point and evaluate f // fullVector<double> xyz(3); for(size_t i = 0; i < nPoint; i++){ xyz(0) = point(i, 0); xyz(1) = point(i, 1); xyz(2) = point(i, 2); eval(i, 0) = f(xyz); } }
double taylorDistanceSq1D(const GradientBasis *gb, const fullMatrix<double> &nodesXYZ, const std::vector<SVector3> &tanCAD) { const int nV = nodesXYZ.size1(); fullMatrix<double> dxyzdX(nV, 3); gb->getGradientsFromNodes(nodesXYZ, &dxyzdX, 0, 0); // const double dx = nodesXYZ(1, 0) - nodesXYZ(0, 0), dy = nodesXYZ(1, 1) - nodesXYZ(0, 1), // dz = nodesXYZ(1, 2) - nodesXYZ(0, 2), h = 0.5*sqrt(dx*dx+dy*dy+dz*dz)/double(nV-1); double distSq = 0.; for (int i=0; i<nV; i++) { SVector3 tanMesh(dxyzdX(i, 0), dxyzdX(i, 1), dxyzdX(i, 2)); const double h = 0.25*tanMesh.normalize(); // Half of "local edge length" SVector3 diff = (dot(tanCAD[i], tanMesh) > 0) ? tanCAD[i] - tanMesh : tanCAD[i] + tanMesh; distSq += h*h*diff.normSq(); } return distSq; }
double l2Norm(const fullMatrix<double>& valA, const fullMatrix<double>& valB){ const size_t nPoint = valA.size1(); const size_t dim = valA.size2(); double norm = 0; double modSquare = 0; for(size_t i = 0; i < nPoint; i++){ modSquare = 0; for(size_t j = 0; j < dim; j++) modSquare += (valA(i, j) - valB(i, j)) * (valA(i, j) - valB(i, j)); norm += modSquare; } return norm = sqrt(norm); }
double taylorDistanceSq2D(const GradientBasis *gb, const fullMatrix<double> &nodesXYZ, const std::vector<SVector3> &normCAD) { const int nV = nodesXYZ.size1(); fullMatrix<double> dxyzdX(nV, 3), dxyzdY(nV, 3); gb->getGradientsFromNodes(nodesXYZ, &dxyzdX, &dxyzdY, 0); double distSq = 0.; for (int i=0; i<nV; i++) { const double nz = dxyzdX(i, 0) * dxyzdY(i, 1) - dxyzdX(i, 1) * dxyzdY(i, 0); const double ny = -dxyzdX(i, 0) * dxyzdY(i, 2) + dxyzdX(i, 2) * dxyzdY(i, 0); const double nx = dxyzdX(i, 1) * dxyzdY(i, 2) - dxyzdX(i, 2) * dxyzdY(i, 1); SVector3 normMesh(nx, ny, nz); const double h = 0.25*sqrt(normMesh.normalize()); // Half of sqrt of "local area", to be adjusted w.r.t. el. type? SVector3 diff = (dot(normCAD[i], normMesh) > 0) ? normCAD[i] - normMesh : normCAD[i] + normMesh; distSq += h*h*diff.normSq(); } return distSq; }
void pyramidalBasis::f(const fullMatrix<double> &coord, fullMatrix<double> &sf) const { const int N = bergot->size(), NPts = coord.size1(); sf.resize(NPts, N); double *fval = new double[N]; for (int iPt=0; iPt<NPts; iPt++) { bergot->f(coord(iPt,0), coord(iPt,1), coord(iPt,2), fval); for (int i = 0; i < N; i++) { sf(iPt,i) = 0.; for (int j = 0; j < N; j++) sf(iPt,i) += coefficients(i,j)*fval[j]; } } delete[] fval; }
void pyramidalBasis::df(const fullMatrix<double> &coord, fullMatrix<double> &dfm) const { const int N = bergot->size(), NPts = coord.size1(); double (*dfv)[3] = new double[N][3]; dfm.resize (N, 3*NPts, false); for (int iPt=0; iPt<NPts; iPt++) { df(coord(iPt,0), coord(iPt,1), coord(iPt,2), dfv); for (int i = 0; i < N; i++) { dfm(i, 3*iPt) = dfv[i][0]; dfm(i, 3*iPt+1) = dfv[i][1]; dfm(i, 3*iPt+2) = dfv[i][2]; } } delete[] dfv; }
void write(bool isScalar, const fullMatrix<double>& l2, string name){ // Stream ofstream stream; string fileName; if(isScalar) fileName = name + "Node.m"; else fileName = name + "Edge.m"; stream.open(fileName.c_str()); // Matrix data const size_t l2Row = l2.size1(); const size_t l2ColMinus = l2.size2() - 1; // Clean Octave stream << "close all;" << endl << "clear all;" << endl << endl; // Mesh (Assuming uniform refinement) stream << "h = [1, "; for(size_t i = 1; i < l2ColMinus; i++) stream << 1 / pow(2, i) << ", "; stream << 1 / pow(2, l2ColMinus) << "];" << endl; // Order (Assuming uniform refinement) stream << "p = [1:" << l2Row << "];" << endl << endl; // Matrix of l2 error (l2[Order][Mesh]) stream << "l2 = ..." << endl << " [..." << endl; for(size_t i = 0; i < l2Row; i++){ stream << " "; for(size_t j = 0; j < l2ColMinus; j++) stream << scientific << showpos << l2(i, j) << " , "; stream << scientific << showpos << l2(i, l2ColMinus) << " ; ..." << endl; } stream << " ];" << endl << endl; // Delta stream << "P = size(p, 2);" << endl << "H = size(h, 2);" << endl << endl << "delta = zeros(P, H - 1);" << endl << endl << "for i = 1:H-1" << endl << " delta(:, i) = ..." << endl << " (log10(l2(:, i + 1)) - log10(l2(:, i))) / ..." << endl << " (log10(1/h(i + 1)) - log10(1/h(i)));" << endl << "end" << endl << endl << "delta" << endl << endl; // Plot stream << "figure;" << endl << "loglog(1./h, l2', '-*');" << endl << "grid;" << endl; // Title stream << "title(" << "'" << name << ": "; if(isScalar) stream << "Nodal"; else stream << "Edge"; stream << "');" << endl << endl; // Axis stream << "xlabel('1/h [-]');" << endl << "ylabel('L2 Error [-]');" << endl; // Close stream.close(); }