Real determinant(const Matrix& m) { #if !defined(QL_NO_UBLAS_SUPPORT) QL_REQUIRE(m.rows() == m.columns(), "matrix is not square"); boost::numeric::ublas::matrix<Real> a(m.rows(), m.columns()); std::copy(m.begin(), m.end(), a.data().begin()); // lu decomposition boost::numeric::ublas::permutation_matrix<Size> pert(m.rows()); /* const Size singular = */ lu_factorize(a, pert); Real retVal = 1.0; for (Size i=0; i < m.rows(); ++i) { if (pert[i] != i) retVal *= -a(i,i); else retVal *= a(i,i); } return retVal; #else QL_FAIL("this version of gcc does not support " "the Boost uBLAS library"); #endif }
Disposable<Matrix> inverse(const Matrix& m) { #if !defined(QL_NO_UBLAS_SUPPORT) QL_REQUIRE(m.rows() == m.columns(), "matrix is not square"); boost::numeric::ublas::matrix<Real> a(m.rows(), m.columns()); std::copy(m.begin(), m.end(), a.data().begin()); boost::numeric::ublas::permutation_matrix<Size> pert(m.rows()); // lu decomposition const Size singular = lu_factorize(a, pert); QL_REQUIRE(singular == 0, "singular matrix given"); boost::numeric::ublas::matrix<Real> inverse = boost::numeric::ublas::identity_matrix<Real>(m.rows()); // backsubstitution boost::numeric::ublas::lu_substitute(a, pert, inverse); Matrix retVal(m.rows(), m.columns()); std::copy(inverse.data().begin(), inverse.data().end(), retVal.begin()); return retVal; #else QL_FAIL("this version of gcc does not support " "the Boost uBLAS library"); #endif }
void TransportGradientPeriodic :: computeTangent(FloatMatrix &k, TimeStep *tStep) { EModelDefaultEquationNumbering fnum; //DofIDEquationNumbering pnum(true, this->grad_ids); EModelDefaultPrescribedEquationNumbering pnum; int nsd = this->domain->giveNumberOfSpatialDimensions(); EngngModel *rve = this->giveDomain()->giveEngngModel(); ///@todo Get this from engineering model std :: unique_ptr< SparseLinearSystemNM > solver( classFactory.createSparseLinSolver( ST_Petsc, this->domain, this->domain->giveEngngModel() ) ); // = rve->giveLinearSolver(); SparseMtrxType stype = solver->giveRecommendedMatrix(true); std :: unique_ptr< SparseMtrx > Kff( classFactory.createSparseMtrx( stype ) ); std :: unique_ptr< SparseMtrx > Kfp( classFactory.createSparseMtrx( stype ) ); std :: unique_ptr< SparseMtrx > Kpp( classFactory.createSparseMtrx( stype ) ); Kff->buildInternalStructure(rve, this->domain->giveNumber(), fnum); int neq = Kff->giveNumberOfRows(); Kfp->buildInternalStructure(rve, this->domain->giveNumber(), fnum, pnum); Kpp->buildInternalStructure(rve, this->domain->giveNumber(), pnum); //Kfp->buildInternalStructure(rve, neq, nsd, {}, {}); //Kpp->buildInternalStructure(rve, nsd, nsd, {}, {}); #if 0 rve->assemble(*Kff, tStep, TangentAssembler(TangentStiffness), fnum, this->domain); rve->assemble(*Kfp, tStep, TangentAssembler(TangentStiffness), fnum, pnum, this->domain); rve->assemble(*Kpp, tStep, TangentAssembler(TangentStiffness), pnum, this->domain); #else auto ma = TangentAssembler(TangentStiffness); IntArray floc, ploc; FloatMatrix mat, R; int nelem = domain->giveNumberOfElements(); #ifdef _OPENMP #pragma omp parallel for shared(Kff, Kfp, Kpp) private(mat, R, floc, ploc) #endif for ( int ielem = 1; ielem <= nelem; ielem++ ) { Element *element = domain->giveElement(ielem); // skip remote elements (these are used as mirrors of remote elements on other domains // when nonlocal constitutive models are used. They introduction is necessary to // allow local averaging on domains without fine grain communication between domains). if ( element->giveParallelMode() == Element_remote || !element->isActivated(tStep) ) { continue; } ma.matrixFromElement(mat, *element, tStep); if ( mat.isNotEmpty() ) { ma.locationFromElement(floc, *element, fnum); ma.locationFromElement(ploc, *element, pnum); ///@todo This rotation matrix is not flexible enough.. it can only work with full size matrices and doesn't allow for flexibility in the matrixassembler. if ( element->giveRotationMatrix(R) ) { mat.rotatedWith(R); } #ifdef _OPENMP #pragma omp critical #endif { Kff->assemble(floc, mat); Kfp->assemble(floc, ploc, mat); Kpp->assemble(ploc, mat); } } } Kff->assembleBegin(); Kfp->assembleBegin(); Kpp->assembleBegin(); Kff->assembleEnd(); Kfp->assembleEnd(); Kpp->assembleEnd(); #endif FloatMatrix grad_pert(nsd, nsd), rhs, sol(neq, nsd); grad_pert.resize(nsd, nsd); grad_pert.beUnitMatrix(); // Workaround since the matrix size is inflexible with custom dof numbering (so far, planned to be fixed). IntArray grad_loc; this->grad->giveLocationArray(this->grad_ids, grad_loc, pnum); FloatMatrix pert(Kpp->giveNumberOfRows(), nsd); pert.assemble(grad_pert, grad_loc, {1,2,3}); //pert.printYourself("pert"); //printf("Kfp = %d x %d\n", Kfp->giveNumberOfRows(), Kfp->giveNumberOfColumns()); //printf("Kff = %d x %d\n", Kff->giveNumberOfRows(), Kff->giveNumberOfColumns()); //printf("Kpp = %d x %d\n", Kpp->giveNumberOfRows(), Kpp->giveNumberOfColumns()); // Compute the solution to each of the pertubation of eps Kfp->times(pert, rhs); //rhs.printYourself("rhs"); // Initial guess (Taylor assumption) helps KSP-iterations for ( auto &n : domain->giveDofManagers() ) { int k1 = n->giveDofWithID( this->dofs(0) )->__giveEquationNumber(); if ( k1 ) { FloatArray *coords = n->giveCoordinates(); for ( int i = 1; i <= nsd; ++i ) { sol.at(k1, i) = -(coords->at(i) - mCenterCoord.at(i)); } } } if ( solver->solve(*Kff, rhs, sol) & NM_NoSuccess ) { OOFEM_ERROR("Failed to solve Kff"); } // Compute the solution to each of the pertubation of eps Kfp->timesT(sol, k); // Assuming symmetry of stiffness matrix // This is probably always zero, but for generality FloatMatrix tmpMat; Kpp->times(pert, tmpMat); k.subtract(tmpMat); k.times( - 1.0 / ( this->domainSize(this->giveDomain(), this->set) + this->domainSize(this->giveDomain(), this->masterSet) )); // Temp workaround on sizing issue mentioned above: FloatMatrix k2 = k; k.beSubMatrixOf(k2, grad_loc, {1,2,3}); }