AnyType MatrixUpWind3::operator()(Stack stack) const { Matrice_Creuse<R> * sparce_mat =GetAny<Matrice_Creuse<R>* >((*emat)(stack)); MatriceMorse<R> * amorse =0; MeshPoint *mp(MeshPointStack(stack)) , mps=*mp; Mesh3 * pTh = GetAny<pmesh3>((*expTh)(stack)); ffassert(pTh); Mesh3 & Th (*pTh); { map< pair<int,int>, R> Aij; KN<double> cc(Th.nv); double infini=DBL_MAX; cc=infini; for (int it=0;it<Th.nt;it++) for (int iv=0;iv<4;iv++) { int i=Th(it,iv); if ( cc[i]==infini) { // if nuset the set mp->setP(&Th,it,iv); cc[i]=GetAny<double>((*expc)(stack)); } } for (int k=0;k<Th.nt;k++) { const Mesh3::Element & K(Th[k]); const Mesh3::Vertex & A(K[0]), &B(K[1]),&C(K[2]),&D(K[3]); R3 Pt(1./4.,1./4.,1./4.); R3 U; MeshPointStack(stack)->set(Th,K(Pt),Pt,K,K.lab); U.x = GetAny< R>( (*expu1)(stack) ) ; U.y = GetAny< R>( (*expu2)(stack) ) ; U.z = GetAny< R>( (*expu3)(stack) ) ; int ii[4] ={ Th(A), Th(B),Th(C),Th(D)};// number of 4 vertex double c[4]={cc[ii[0]],cc[ii[1]],cc[ii[2]],cc[ii[3]]}; double a[4][4]; if (Marco(K,U,c,a) ) { for (int i=0;i<4;i++) for (int j=0;j<4;j++) if (fabs(a[i][j]) >= 1e-30) Aij[make_pair(ii[i],ii[j])]+=a[i][j]; } } amorse= new MatriceMorse<R>(Th.nv,Th.nv,Aij,false); } sparce_mat->Uh=UniqueffId(); sparce_mat->Vh=UniqueffId(); sparce_mat->A.master(amorse); sparce_mat->typemat=(amorse->n == amorse->m) ? TypeSolveMat(TypeSolveMat::GMRES) : TypeSolveMat(TypeSolveMat::NONESQUARE); // none square matrice (morse) *mp=mps; if(verbosity>3) { cout << " End Build MatrixUpWind : " << endl;} return sparce_mat; }
AnyType MatrixUpWind0::operator()(Stack stack) const { Matrice_Creuse<R> * sparce_mat =GetAny<Matrice_Creuse<R>* >((*emat)(stack)); MatriceMorse<R> * amorse =0; MeshPoint *mp(MeshPointStack(stack)) , mps=*mp; Mesh * pTh = GetAny<pmesh>((*expTh)(stack)); ffassert(pTh); Mesh & Th (*pTh); { map< pair<int,int>, R> Aij; KN<double> cc(Th.nv); double infini=DBL_MAX; cc=infini; for (int it=0;it<Th.nt;it++) for (int iv=0;iv<3;iv++) { int i=Th(it,iv); if ( cc[i]==infini) { // if nuset the set mp->setP(&Th,it,iv); cc[i]=GetAny<double>((*expc)(stack)); } } for (int k=0;k<Th.nt;k++) { const Triangle & K(Th[k]); const Vertex & A(K[0]), &B(K[1]),&C(K[2]); R2 Pt(1./3.,1./3.); R u[2]; MeshPointStack(stack)->set(Th,K(Pt),Pt,K,K.lab); u[0] = GetAny< R>( (*expu1)(stack) ) ; u[1] = GetAny< R>( (*expu2)(stack) ) ; int ii[3] ={ Th(A), Th(B),Th(C)}; double q[3][2]= { { A.x,A.y} ,{B.x,B.y},{C.x,C.y} } ; // coordinates of 3 vertices (input) double c[3]={cc[ii[0]],cc[ii[1]],cc[ii[2]]}; double a[3][3]; if (gladys(q,u,c,a) ) { for (int i=0;i<3;i++) for (int j=0;j<3;j++) if (fabs(a[i][j]) >= 1e-30) Aij[make_pair(ii[i],ii[j])]+=a[i][j]; } } amorse= new MatriceMorse<R>(Th.nv,Th.nv,Aij,false); } sparce_mat->Uh=UniqueffId(); sparce_mat->Vh=UniqueffId(); sparce_mat->A.master(amorse); sparce_mat->typemat=(amorse->n == amorse->m) ? TypeSolveMat(TypeSolveMat::GMRES) : TypeSolveMat(TypeSolveMat::NONESQUARE); // none square matrice (morse) *mp=mps; if(verbosity>3) { cout << " End Build MatrixUpWind : " << endl;} return sparce_mat; }
void Mesh2::GRead(FILE * ff) { PlotStream f(ff); string s; f >> s; ffassert( s== Gsbegin); f >> nv >> nt >> nbe; if(verbosity>1) cout << " GRead : nv " << nv << " " << nt << " " << nbe << endl; this->vertices = new Vertex[nv]; this->elements = new Element [nt]; this->borderelements = new BorderElement[nbe]; for (int k=0; k<nv; k++) { Vertex & P = this->vertices[k]; f >> P.x >>P.y >> P.lab ; } mes=0.; mesb=0.; if(nt != 0) { for (int k=0; k<nt; k++) { int i[4],lab; Element & K(this->elements[k]); f >> i[0] >> i[1] >> i[2] >> lab; K.set(this->vertices,i,lab); mes += K.mesure(); } } for (int k=0; k<nbe; k++) { int i[4],lab; BorderElement & K(this->borderelements[k]); f >> i[0] >> i[1] >> lab; K.set(this->vertices,i,lab); mesb += K.mesure(); } f >> s; ffassert( s== Gsend); }
void exchange_restriction(Type* const& pA, KN<K>* pin, KN<K>* pout, MatriceMorse<double>* mR) { if(pA->_exchange && !pA->_exchange[1]) { ffassert((!mR && pA->_exchange[0]->getDof() == pout->n) || (mR && mR->n == pin->n && mR->m == pout->n)); PETSc::changeNumbering_func(pA, pin, pout, false); PETSc::changeNumbering_func(pA, pin, pout, true); pout->resize(pA->_exchange[0]->getDof()); *pout = K(); if(mR) { for(int i = 0; i < mR->n; ++i) { for(int j = mR->lg[i]; j < mR->lg[i + 1]; ++j) pout->operator[](mR->cl[j]) += mR->a[j] * pin->operator[](i); } } exchange_dispatched(pA->_exchange[0], pout, false); } }
AnyType exchangeInOut_Op<Type, K>::operator()(Stack stack) const { Type* pA = GetAny<Type*>((*A)(stack)); KN<K>* pin = GetAny<KN<K>*>((*in)(stack)); KN<K>* pout = GetAny<KN<K>*>((*out)(stack)); bool scaled = nargs[0] && GetAny<bool>((*nargs[0])(stack)); Matrice_Creuse<double>* pR = nargs[1] ? GetAny<Matrice_Creuse<double>*>((*nargs[1])(stack)) : nullptr; MatriceMorse<double>* mR = pR ? static_cast<MatriceMorse<double>*>(&(*pR->A)) : nullptr; if(pR) { ffassert(!scaled); exchange_restriction(pA, pin, pout, mR); } else if(pin->n == pout->n) { *pout = *pin; exchange(pA, pout, scaled); } return 0L; }
void TypeOfFE_PkEdge::FB(const bool * whatd,const Mesh & ,const Triangle & K,const R2 & P,RNMK_ & val) const { R2 A(K[0]), B(K[1]),C(K[2]); R l0=1-P.x-P.y,l1=P.x,l2=P.y; R L[3]={l0,l1,l2}; assert( val.N()>=ndf); assert(val.M()==1); int ee=0; if (L[0] <= min(L[1],L[2]) ) ee=0; // arete else if (L[1] <= min(L[0],L[2]) ) ee=1; else ee=2; int e3=ee*npe; double s=1.-L[ee]; R xe = L[VerticesOfTriangularEdge[ee][0]]/s;// go from 0 to 1 on edge if(K.EdgeOrientation(ee) <0.) xe = 1-xe; //cout << P << " ee = " << ee << " xe " << xe << " " << L[ee]<< " s=" <<s << " orient: " << K.EdgeOrientation(ee) <<endl; assert(s); val=0; if (whatd[op_id]) { RN_ f0(val('.',0,op_id)); for (int l=0;l<npe;l++) { int df= e3+l; R f=1.; for (int i=0;i<npe;++i) if(i != l) f *= (xe-X[i])/(X[l]-X[i]); f0[df] = f; } //cout << " f0 = " << f0 << " X= "<< X << endl; } if( whatd[op_dx] || whatd[op_dy] || whatd[op_dxx] || whatd[op_dyy] || whatd[op_dxy]) { cerr << " TO DO ??? FH " << endl; ffassert(0); } }
void addMatMul (const KN_<R> &x, KN_<R> &Ax) const { ffassert(x.N() == Ax.N()); Ax += (const MatriceMorse<R> &)(*this) * x; }
void Solver (const MatriceMorse<R> &AA, KN_<R> &x, const KN_<R> &b) const { SuperMatrix B, X; SuperLUStat_t stat; int info = 0, lwork = 0; double ferr[1], berr[1]; double rpg, rcond; double *xx; B.Store = 0; X.Store = 0; ffassert(&x[0] != &b[0]); epsr = (eps < 0) ? (epsr > 0 ? -epsr : -eps) : eps; Dtype_t R_SLU = SuperLUDriver<R>::R_SLU_T(); { void *work = 0; int nrhs = 1; KN_2Ptr<R> xx(x), bb(b); // cout << " xx #### " << xx.c.N() << " "<< xx.ca.N() << " " << xx.ca.step << endl; // cout << " bb #### " << bb.c.N() << " "<< bb.ca.N() << " " << bb.ca.step <<endl; // FFCS - "this->" required by g++ 4.7 this->Create_Dense_Matrix(&B, m, 1, bb, m, SLU_DN, R_SLU, SLU_GE); this->Create_Dense_Matrix(&X, m, 1, xx, m, SLU_DN, R_SLU, SLU_GE); B.ncol = nrhs; /* Set the number of right-hand side */ /* Initialize the statistics variables. */ StatInit(&stat); SuperLUDriver<R>::gssvx(&options, &A, perm_c, perm_r, etree, equed, RR, CC, &L, &U, work, lwork, &B, &X, &rpg, &rcond, ferr, berr, &Glu, &mem_usage, &stat, &info); if (verbosity > 2) { printf("Triangular solve: dgssvx() returns info %d\n", info); } } if (verbosity > 3) { if (info == 0 || info == n + 1) { /* This is how you could access the solution matrix. */ // R *sol = (R *)((DNformat *)X.Store)->nzval; if (options.IterRefine) { int i = 0; printf("Iterative Refinement:\n"); printf("%8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); printf("%8d%8d%16e%16e\n", i + 1, stat.RefineSteps, ferr[0], berr[0]); } fflush(stdout); } else if (info > 0 && lwork == -1) { printf("** Estimated memory: %d bytes\n", info - n); } } // cout << " x min max " << x.min() << " " <<x.max() << endl; // cout << "=========================================" << endl; if (B.Store) {Destroy_SuperMatrix_Store(&B);} if (X.Store) {Destroy_SuperMatrix_Store(&X);} }
AnyType removeDOF_Op<T>::operator()(Stack stack) const { Matrice_Creuse<T>* pA = GetAny<Matrice_Creuse<T>* >((*A)(stack)); Matrice_Creuse<T>* pR = GetAny<Matrice_Creuse<T>* >((*R)(stack)); KN<T>* pX = GetAny<KN<T>* >((*x)(stack)); KN<T>* pOut = GetAny<KN<T>* >((*out)(stack)); ffassert(pA && pR && pX && pOut); pA->Uh = pR->Uh; pA->Vh = pR->Vh; MatriceMorse<T> *mA = static_cast<MatriceMorse<T>*>(&(*pA->A)); MatriceMorse<T> *mR = static_cast<MatriceMorse<T>*>(&(*pR->A)); bool symmetrize = nargs[0] ? GetAny<bool>((*nargs[0])(stack)) : false; KN<long>* condensation = nargs[1] ? GetAny<KN<long>* >((*nargs[1])(stack)) : (KN<long>*) 0; unsigned int n = condensation ? condensation->n : mR->nbcoef; int* lg = new int[n + 1]; int* cl; T* val; T* b; if(pOut->n == 0) { b = new T[n]; pOut->set(b, n); } std::vector<signed int> tmpVec; if(!condensation) { tmpVec.resize(mA->n); for(unsigned int i = 0; i < n; ++i) tmpVec[mR->cl[i]] = i + 1; if(!mA->symetrique) { std::vector<std::pair<int, T> > tmp; tmp.reserve(mA->nbcoef); lg[0] = 0; for(unsigned int i = 0; i < n; ++i) { for(unsigned int j = mA->lg[mR->cl[i]]; j < mA->lg[mR->cl[i] + 1]; ++j) { unsigned int col = tmpVec[mA->cl[j]]; if(col != 0 && abs(mA->a[j]) > EPS) { if(symmetrize) { if(col - 1 <= i) tmp.emplace_back(col - 1, mA->a[j]); } else tmp.emplace_back(col - 1, mA->a[j]); } } std::sort(tmp.begin() + lg[i], tmp.end(), [](const std::pair<unsigned int, T>& lhs, const std::pair<unsigned int, T>& rhs) { return lhs.first < rhs.first; }); *(*pOut + i) = *(*pX + mR->cl[i]); lg[i + 1] = tmp.size(); } mA->nbcoef = tmp.size(); if(symmetrize) mA->symetrique = true; else mA->symetrique = false; cl = new int[tmp.size()]; val = new T[tmp.size()]; for(unsigned int i = 0; i < tmp.size(); ++i) { cl[i] = tmp[i].first; val[i] = tmp[i].second; } } else { std::vector<std::vector<std::pair<unsigned int, T> > > tmp(n); for(unsigned int i = 0; i < n; ++i) tmp[i].reserve(mA->lg[mR->cl[i] + 1] - mA->lg[mR->cl[i]]); unsigned int nnz = 0; for(unsigned int i = 0; i < n; ++i) { for(unsigned int j = mA->lg[mR->cl[i]]; j < mA->lg[mR->cl[i] + 1]; ++j) { unsigned int col = tmpVec[mA->cl[j]]; if(col != 0 && abs(mA->a[j]) > EPS) { if(i < col - 1) tmp[col - 1].emplace_back(i, mA->a[j]); else tmp[i].emplace_back(col - 1, mA->a[j]); ++nnz; } } *(*pOut + i) = *(*pX + mR->cl[i]); } mA->nbcoef = nnz; cl = new int[nnz]; val = new T[nnz]; nnz = 0; lg[0] = 0; for(unsigned int i = 0; i < n; ++i) { std::sort(tmp[i].begin(), tmp[i].end(), [](const std::pair<unsigned int, T>& lhs, const std::pair<unsigned int, T>& rhs) { return lhs.first < rhs.first; }); for(typename std::vector<std::pair<unsigned int, T> >::const_iterator it = tmp[i].begin(); it != tmp[i].end(); ++it) { cl[nnz] = it->first; val[nnz++] = it->second; } lg[i + 1] = nnz; } } delete [] mA->cl; delete [] mA->lg; delete [] mA->a; mA->n = n; mA->m = n; mA->N = n; mA->M = n; mA->lg = lg; mA->cl = cl; mA->a = val; } else { tmpVec.reserve(mA->n); unsigned int i = 0, j = 1; for(unsigned int k = 0; k < mA->n; ++k) { if(k == *(*condensation + i)) { ++i; tmpVec.emplace_back(i); } else { tmpVec.emplace_back(-j); ++j; } } // if(!mA->symetrique) { std::vector<std::pair<int, T> > tmpInterior; std::vector<std::pair<int, T> > tmpBoundary; std::vector<std::pair<int, T> > tmpInteraction; tmpInterior.reserve(mA->nbcoef); tmpBoundary.reserve(mA->nbcoef); tmpInteraction.reserve(mA->nbcoef); lg[0] = 0; for(unsigned int i = 0; i < mA->n; ++i) { int row = tmpVec[i]; if(row < 0) { for(unsigned int j = mA->lg[i]; j < mA->lg[i + 1]; ++j) { int col = tmpVec[mA->cl[j]]; if(col < 0) tmpInterior.emplace_back(-col - 1, mA->a[j]); else tmpInteraction.emplace_back(col - 1, mA->a[j]); } } else { for(unsigned int j = mA->lg[i]; j < mA->lg[i + 1]; ++j) { int col = tmpVec[mA->cl[j]]; if(col > 0) tmpBoundary.emplace_back(col - 1, mA->a[j]); } // std::sort(tmp.begin() + lg[i], tmp.end()); *(*pOut + i) = *(*pX + *(*condensation + i)); lg[i + 1] = tmpBoundary.size(); } } cl = new int[tmpBoundary.size()]; val = new T[tmpBoundary.size()]; for(unsigned int i = 0; i < tmpBoundary.size(); ++i) { cl[i] = tmpBoundary[i].first; val[i] = tmpBoundary[i].second; } // } MatriceMorse<T>* m = new MatriceMorse<T>(n, n, tmpBoundary.size(), mA->symetrique, val, lg, cl, true); pR->typemat = TypeSolveMat(TypeSolveMat::GMRES); pR->A.master(m); m->dummy = false; } return 0L; }