ordinal_type
Stokhos::MonomialGramSchmidtPCEBasis<ordinal_type, value_type>::
buildReducedBasis(
  ordinal_type max_p, 
  value_type threshold, 
  const Teuchos::SerialDenseMatrix<ordinal_type,value_type>& A, 
  const Teuchos::SerialDenseMatrix<ordinal_type,value_type>& F,
  const Teuchos::Array<value_type>& weights, 
  Teuchos::Array< Stokhos::MultiIndex<ordinal_type> >& terms_,
  Teuchos::Array<ordinal_type>& num_terms_,
  Teuchos::SerialDenseMatrix<ordinal_type,value_type>& Qp_, 
  Teuchos::SerialDenseMatrix<ordinal_type,value_type>& Q_)
{
  // Compute basis terms -- 2-D array giving powers for each linear index
  ordinal_type max_sz;
  CPBUtils::compute_terms(max_p, this->d, max_sz, terms_, num_terms_);

  // Compute B matrix -- monomials in F
  // for i=0,...,nqp-1
  //   for j=0,...,sz-1
  //      B(i,j) = F(i,1)^terms_[j][1] * ... * F(i,d)^terms_[j][d]
  // where sz is the total size of a basis up to order p and terms_[j] 
  // is an array of powers for each term in the total-order basis
  ordinal_type nqp = weights.size();
  SDM B(nqp, max_sz);
  for (ordinal_type i=0; i<nqp; i++) {
    for (ordinal_type j=0; j<max_sz; j++) {
      B(i,j) = 1.0;
      for (ordinal_type k=0; k<this->d; k++)
	B(i,j) *= std::pow(F(i,k), terms_[j][k]);
    }
  }

  // Rescale columns of B to have unit norm
  for (ordinal_type j=0; j<max_sz; j++) {
    value_type nrm = 0.0;
    for (ordinal_type i=0; i<nqp; i++)
      nrm += B(i,j)*B(i,j)*weights[i];
    nrm = std::sqrt(nrm);
    for (ordinal_type i=0; i<nqp; i++)
      B(i,j) /= nrm;
  }

  // Compute our new basis -- each column of Q is the new basis evaluated
  // at the original quadrature points.  Constraint pivoting so first d+1
  // columns and included in Q.
  SDM R;
  Teuchos::Array<ordinal_type> piv(max_sz);
  for (int i=0; i<this->d+1; i++)
    piv[i] = 1;
  typedef Stokhos::OrthogonalizationFactory<ordinal_type,value_type> SOF;
   ordinal_type sz_ = SOF::createOrthogonalBasis(
    this->orthogonalization_method, threshold, this->verbose, B, weights, 
    Q_, R, piv);

  // Compute Qp = A^T*W*Q
  SDM tmp(nqp, sz_);
  Qp_.reshape(this->pce_sz, sz_);
  for (ordinal_type i=0; i<nqp; i++)
    for (ordinal_type j=0; j<sz_; j++)
      tmp(i,j) = Q_(i,j)*weights[i];
  ordinal_type ret = 
    Qp_.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, A, tmp, 0.0);
  TEUCHOS_ASSERT(ret == 0);

  // It isn't clear that Qp is orthogonal, but if you derive the projection
  // matrix from the original space to the reduced, you end up with 
  // Q^T*W*A = Qp^T

  return sz_;
}
Exemplo n.º 2
0
ordinal_type
Stokhos::MonomialProjGramSchmidtPCEBasis<ordinal_type, value_type>::
buildReducedBasis(
  ordinal_type max_p, 
  value_type threshold,
  const Teuchos::SerialDenseMatrix<ordinal_type,value_type>& A, 
  const Teuchos::SerialDenseMatrix<ordinal_type,value_type>& F,
  const Teuchos::Array<value_type>& weights, 
  Teuchos::Array< Stokhos::MultiIndex<ordinal_type> >& terms_,
  Teuchos::Array<ordinal_type>& num_terms_,
  Teuchos::SerialDenseMatrix<ordinal_type,value_type>& Qp_, 
  Teuchos::SerialDenseMatrix<ordinal_type,value_type>& Q_)
{
  // Compute basis terms -- 2-D array giving powers for each linear index
  ordinal_type max_sz;
  CPBUtils::compute_terms(max_p, this->d, max_sz, terms_, num_terms_);

  // Compute B matrix -- monomials in F
  // for i=0,...,nqp-1
  //   for j=0,...,sz-1
  //      B(i,j) = F(i,1)^terms_[j][1] * ... * F(i,d)^terms_[j][d]
  // where sz is the total size of a basis up to order p and terms_[j] 
  // is an array of powers for each term in the total-order basis
  ordinal_type nqp = weights.size();
  SDM B(nqp, max_sz);
  for (ordinal_type i=0; i<nqp; i++) {
    for (ordinal_type j=0; j<max_sz; j++) {
      B(i,j) = 1.0;
      for (ordinal_type k=0; k<this->d; k++)
	B(i,j) *= std::pow(F(i,k), terms_[j][k]);
    }
  }

  // Project B into original basis -- should use SPAM for this
  SDM Bp(this->pce_sz, max_sz);
  const Teuchos::Array<value_type>& basis_norms = 
    this->pce_basis->norm_squared();
  for (ordinal_type i=0; i<this->pce_sz; i++) {
    for (ordinal_type j=0; j<max_sz; j++) {
      Bp(i,j) = 0.0;
      for (ordinal_type k=0; k<nqp; k++)
	Bp(i,j) += weights[k]*B(k,j)*A(k,i);
      Bp(i,j) /= basis_norms[i];
    }
  }

  // Rescale columns of Bp to have unit norm
  for (ordinal_type j=0; j<max_sz; j++) {
    value_type nrm = 0.0;
    for (ordinal_type i=0; i<this->pce_sz; i++)
      nrm += Bp(i,j)*Bp(i,j)*basis_norms[i];
    nrm = std::sqrt(nrm);
    for (ordinal_type i=0; i<this->pce_sz; i++)
      Bp(i,j) /= nrm;
  }

  // Compute our new basis -- each column of Qp is the coefficients of the
  // new basis in the original basis.  Constraint pivoting so first d+1
  // columns and included in Qp.
  Teuchos::Array<value_type> w(this->pce_sz, 1.0);
  SDM R;
  Teuchos::Array<ordinal_type> piv(max_sz);
  for (int i=0; i<this->d+1; i++)
    piv[i] = 1;
  typedef Stokhos::OrthogonalizationFactory<ordinal_type,value_type> SOF;
  ordinal_type sz_ = SOF::createOrthogonalBasis(
    this->orthogonalization_method, threshold, this->verbose, Bp, w, 
    Qp_, R, piv);

  // Evaluate new basis at original quadrature points
  Q_.reshape(nqp, sz_);
  ordinal_type ret = 
    Q_.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, A, Qp_, 0.0);
  TEUCHOS_ASSERT(ret == 0);

  return sz_;
}
Exemplo n.º 3
0
//GMRES  
int gmres(const  Teuchos::SerialDenseMatrix<int, double> &  A, Teuchos::SerialDenseMatrix<int,double>   X,const Teuchos::SerialDenseMatrix<int,double> &   B, int max_iter, double tolerance)

{
  int n; 
  int k;
  double resid;
  k=1;
  n=A.numRows();
  std::cout << "A= " << A << std::endl;
  std::cout << "B= " << B << std::endl;
  //Teuchos::SerialDenseMatrix<int, double> Ax(n,1);
  //Ax.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,1.0, A, X, 0.0);

  Teuchos::SerialDenseMatrix<int, double> r0(B);
  //r0-=Ax;
    
  resid=r0.normFrobenius();
  std::cout << "resid= " << resid << std::endl;
  //define vector v=r/norm(r) where r=b-Ax
  
  r0.scale(1/resid);
  
  Teuchos::SerialDenseMatrix<int, double> h(1,1);

  //Matrix of orthog basis vectors V
  Teuchos::SerialDenseMatrix<int, double> V(n,1);
  
   //Set v=r0/norm(r0) to be 1st col of V
   for (int i=0; i<n; i++){
        V(i,0)=r0(i,0);
       }
   //right hand side
   Teuchos::SerialDenseMatrix<int, double> bb(1,1);
   bb(0,0)=resid;
   Teuchos::SerialDenseMatrix<int, double> w(n,1);
   Teuchos::SerialDenseMatrix<int, double> c;
   Teuchos::SerialDenseMatrix<int, double> s;
  
   while (resid > tolerance && k < max_iter){
    
    std::cout << "k = " << k << std::endl;
    h.reshape(k+1,k);
    //Arnoldi iteration(Gram-Schmidt )
    V.reshape(n,k+1);    
    //set vk to be kth col of V
    Teuchos::SerialDenseMatrix<int, double> vk(Teuchos::Copy, V, n,1,0,k-1);
    
    w.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, A, vk, 0.0);  
    Teuchos::SerialDenseMatrix<int, double> vi(n,1);
    Teuchos::SerialDenseMatrix<int, double> ip(1,1);
    for (int i=0; i<k; i++){
       //set vi to be ith col of V
       Teuchos::SerialDenseMatrix<int, double> vi(Teuchos::Copy, V, n,1,0,i);    
       //Calculate inner product
       ip.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, vi, w, 0.0);
       h(i,k-1)= ip(0,0);
       //scale vi by h(i,k-1)
       vi.scale(ip(0,0));     
       w-=vi;
       }         
    h(k,k-1)=w.normFrobenius();     
 
    w.scale(1.0/w.normFrobenius());   
    //add column vk+1=w to V
    for (int i=0; i<n; i++){
          V(i,k)=w(i,0);
         } 
    
   //Solve upper hessenberg least squares problem via Givens rotations
   //Compute previous Givens rotations
    for (int i=0; i<k-1; i++){
     //  double hi=h(i,k-1);
     //  double hi1=h(i+1,k-1);

     // h(i,k-1)=c(i,0)*h(i,k-1)+s(i,0)*h(i+1,k-1);
     // h(i+1,k-1)=-1*s(i,0)*h(i,k-1)+c(i,0)*h(i+1,k-1);
      // h(i,k-1)=c(i,0)*hi+s(i,0)*hi1;
      // h(i+1,k-1)=-1*s(i,0)*hi+c(i,0)*hi1;   
     
     double q=c(i,0)*h(i,k-1)+s(i,0)*h(i+1,k-1);
     h(i+1,k-1)=-1*s(i,0)*h(i,k-1)+c(i,0)*h(i+1,k-1);
     h(i,k-1)=q;




     }  
     //Compute next Givens rotations
     c.reshape(k,1);
     s.reshape(k,1); 
     bb.reshape(k+1,1);
     double l = sqrt(h(k-1,k-1)*h(k-1,k-1)+h(k,k-1)*h(k,k-1));
     c(k-1,0)=h(k-1,k-1)/l;
     s(k-1,0)=h(k,k-1)/l;
     
     std::cout << "c  "  <<  c(k-1,0)<<std::endl;
     std::cout << "s "  <<  s(k-1,0)<<std::endl;

    
     // Givens rotation on h and bb
       
   //  h(k-1,k-1)=l;
     
  //  h(k,k-1)=0;
       double hk=h(k,k-1);
       double hk1=h(k-1,k-1);

      h(k-1,k-1)=c(k-1,0)*hk1+s(k-1,0)*hk;
      h(k,k-1)=-1*s(k-1,0)*hk1+c(k-1,0)*hk;

     std::cout << "l = " << l <<std::endl;
     std::cout << "h(k-1,k-1) = should be l  " << h(k-1,k-1) <<std::endl;
     std::cout << "h(k,k-1) = should be 0  " << h(k,k-1) <<std::endl;
     bb(k,0)=-1*s(k-1,0)*bb(k-1,0); 
     bb(k-1,0)=c(k-1,0)*bb(k-1,0);
     
   
    //Determine residual    
     resid =fabs(bb(k,0));
      
     std::cout << "resid = " << resid <<std::endl;
     k++;
  } 
  
  //Extract upper triangular square matrix
   bb.reshape(h.numRows()-1 ,1);
   
   //Solve linear system
   int info;
   std::cout  << "bb pre solve = " << bb << std::endl;
   std::cout << "h= " << h << std::endl;
   Teuchos::LAPACK<int, double> lapack;
   lapack.TRTRS('U', 'N', 'N', h.numRows()-1, 1, h.values(), h.stride(), bb.values(), bb.stride(),&info); 

   V.reshape(n,k-1);
   
   std::cout  << "V= " << V << std::endl;
   std::cout  << "y= " << bb << std::endl;
   X.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, V, bb, 1.0);
   std::cout << "X=  " << X << std::endl;

  


   //Check V is orthogoanl
  // Teuchos::SerialDenseMatrix<int, double> vtv(V);
  // vtv.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, V, V, 0.0);
  // std::cout << "Vtv" << vtv << std::endl;

return 0;
}
Exemplo n.º 4
0
int
Stokhos::GMRESDivisionExpansionStrategy<ordinal_type,value_type,node_type>::
GMRES(const Teuchos::SerialDenseMatrix<int, double> &  A, Teuchos::SerialDenseMatrix<int,double> &  X, const Teuchos::SerialDenseMatrix<int,double> &   B, int max_iter, double tolerance, int prec_iter, int order, int dim, int PrecNum, const Teuchos::SerialDenseMatrix<int, double> & M, int diag)
{
    int n = A.numRows();
    int k = 1;
    double resid;
    Teuchos::SerialDenseMatrix<int, double> P(n,n);
    Teuchos::SerialDenseMatrix<int, double> Ax(n,1);
    Ax.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,1.0, A, X, 0.0);
    Teuchos::SerialDenseMatrix<int, double> r0(B);
    r0-=Ax;
    resid=r0.normFrobenius();
//define vector v=r/norm(r) where r=b-Ax
    Teuchos::SerialDenseMatrix<int, double> v(n,1);
    r0.scale(1/resid);
    Teuchos::SerialDenseMatrix<int, double> h(1,1);
//Matrix of orthog basis vectors V
    Teuchos::SerialDenseMatrix<int, double> V(n,1);
//Set v=r0/norm(r0) to be 1st col of V
    for (int i=0; i<n; i++) {
        V(i,0)=r0(i,0);
    }
    //right hand side
    Teuchos::SerialDenseMatrix<int, double> bb(1,1);
    bb(0,0)=resid;
    Teuchos::SerialDenseMatrix<int, double> w(n,1);
    Teuchos::SerialDenseMatrix<int, double> c;
    Teuchos::SerialDenseMatrix<int, double> s;
    while (resid > tolerance && k < max_iter) {
        h.reshape(k+1,k);
        //Arnoldi iteration(Gram-Schmidt )
        V.reshape(n,k+1);
        //set vk to be kth col of V
        Teuchos::SerialDenseMatrix<int, double> vk(Teuchos::Copy, V, n,1,0,k-1);
        //Preconditioning step: solve Mz=vk
        Teuchos::SerialDenseMatrix<int, double> z(vk);
        if (PrecNum == 1) {
            Stokhos::DiagPreconditioner precond(M);
            precond.ApplyInverse(vk,z,prec_iter);
        }
        else if (PrecNum == 2) {
            Stokhos::JacobiPreconditioner precond(M);
            precond.ApplyInverse(vk,z,2);
        }
        else if (PrecNum == 3) {
            Stokhos::GSPreconditioner precond(M,1);
            precond.ApplyInverse(vk,z,1);
        }
        else if (PrecNum == 4) {
            Stokhos::SchurPreconditioner precond(M, order, dim, diag);
            precond.ApplyInverse(vk,z,prec_iter);
        }

        w.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1, A, z, 0.0);
        Teuchos::SerialDenseMatrix<int, double> vi(n,1);
        Teuchos::SerialDenseMatrix<int, double> ip(1,1);
        for (int i=0; i<k; i++) {
            //set vi to be ith col of V
            Teuchos::SerialDenseMatrix<int, double> vi(Teuchos::Copy, V, n,1,0,i);
            //Calculate inner product
            ip.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, vi, w, 0.0);
            h(i,k-1)= ip(0,0);
            //scale vi by h(i,k-1)
            vi.scale(ip(0,0));
            w-=vi;
        }
        h(k,k-1)=w.normFrobenius();
        w.scale(1.0/h(k,k-1));
        //add column vk+1=w to V
        for (int i=0; i<n; i++) {
            V(i,k)=w(i,0);
        }
        //Solve upper hessenberg least squares problem via Givens rotations
        //Compute previous Givens rotations
        for (int i=0; i<k-1; i++) {
            double q=c(i,0)*h(i,k-1)+s(i,0)*h(i+1,k-1);
            h(i+1,k-1)=-1*s(i,0)*h(i,k-1)+c(i,0)*h(i+1,k-1);
            h(i,k-1)=q;

        }
        //Compute next Givens rotations
        c.reshape(k,1);
        s.reshape(k,1);
        bb.reshape(k+1,1);
        double l = sqrt(h(k-1,k-1)*h(k-1,k-1)+h(k,k-1)*h(k,k-1));
        c(k-1,0)=h(k-1,k-1)/l;
        s(k-1,0)=h(k,k-1)/l;

        // Givens rotation on h and bb
        h(k-1,k-1)=l;
        h(k,k-1)=0;

        bb(k,0)=-s(k-1,0)*bb(k-1,0);
        bb(k-1,0)=c(k-1,0)*bb(k-1,0);

        //Determine residual
        resid = fabs(bb(k,0));
        k++;
    }
    //Extract upper triangular square matrix
    bb.reshape(h.numRows()-1 ,1);
    //Solve linear system
    int info;
    Teuchos::LAPACK<int, double> lapack;
    lapack.TRTRS('U', 'N', 'N', h.numRows()-1, 1, h.values(), h.stride(), bb.values(), bb.stride(),&info);
    Teuchos::SerialDenseMatrix<int, double> ans(X);
    V.reshape(n,k-1);
    ans.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, V, bb, 0.0);
    if (PrecNum == 1) {
        Stokhos::DiagPreconditioner precond(M);
        precond.ApplyInverse(ans,ans,prec_iter);
    }
    else if (PrecNum == 2) {
        Stokhos::JacobiPreconditioner precond(M);
        precond.ApplyInverse(ans,ans,2);
    }
    else if (PrecNum == 3) {
        Stokhos::GSPreconditioner precond(M,1);
        precond.ApplyInverse(ans,ans,1);
    }
    else if (PrecNum == 4) {
        Stokhos::SchurPreconditioner precond(M, order, dim, diag);
        precond.ApplyInverse(ans,ans,prec_iter);
    }
    X+=ans;

    std::cout << "iteration count=  " << k-1 << std::endl;


    return 0;
}
Exemplo n.º 5
0
//Mean-Based Preconditioned GMRES  
int pregmres(const Teuchos::SerialDenseMatrix<int, double> &  A, const Teuchos::SerialDenseMatrix<int,double> &  X,const Teuchos::SerialDenseMatrix<int,double> &   B, int max_iter, double tolerance)
{
  int n; 
  int k;
  double resid;
  k=1;
  n=A.numRows();
  std::cout << A << std::endl;
  Teuchos::SerialDenseMatrix<int, double> D(n,1);

  //Get diagonal entries of A 
  for (int i=0; i<n; i++){
    D(i,0)=A(i,i);
  }
  

  Teuchos::SerialDenseMatrix<int, double> Ax(n,1);
  Ax.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,1.0, A, X, 0.0);

  Teuchos::SerialDenseMatrix<int, double> r0(B);
  r0-=Ax;
  
  resid=r0.normFrobenius();
  
  //define vector v=r/norm(r) where r=b-Ax
  Teuchos::SerialDenseMatrix<int, double> v(n,1);
  r0.scale(1/resid);
  
  Teuchos::SerialDenseMatrix<int, double> h(1,1);

  //Matrix of orthog basis vectors V
  Teuchos::SerialDenseMatrix<int, double> V(n,1);
  
   //Set v=r0/norm(r0) to be 1st col of V
   for (int i=0; i<n; i++){
        V(i,0)=r0(i,0);
       }
   //right hand side
   Teuchos::SerialDenseMatrix<int, double> bb(1,1);
   bb(0,0)=resid;
   Teuchos::SerialDenseMatrix<int, double> w(n,1);
   Teuchos::SerialDenseMatrix<int, double> c;
   Teuchos::SerialDenseMatrix<int, double> s;
       
   while (resid > tolerance && k < max_iter){
    
    std::cout << "k = " << k << std::endl;
    h.reshape(k+1,k);
    //Arnoldi iteration(Gram-Schmidt )
    V.reshape(n,k+1);    
    //set vk to be kth col of V
    Teuchos::SerialDenseMatrix<int, double> vk(Teuchos::Copy, V, n,1,0,k-1);
    //Preconditioning step w=AMj(-1)vj
    w.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1/D(k-1,0), A, vk, 0.0);  

   
    
    Teuchos::SerialDenseMatrix<int, double> vi(n,1);
    Teuchos::SerialDenseMatrix<int, double> ip(1,1);
    for (int i=0; i<k; i++){
       //set vi to be ith col of V
       Teuchos::SerialDenseMatrix<int, double> vi(Teuchos::Copy, V, n,1,0,i);    
       //Calculate inner product
       ip.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, vi, w, 0.0);
       h(i,k-1)= ip(0,0);
       //scale vi by h(i,k-1)
       vi.scale(ip(0,0));     
       w-=vi;
       }         
    h(k,k-1)=w.normFrobenius();     
    
    w.scale(1.0/w.normFrobenius());   
    //add column vk+1=w to V
    for (int i=0; i<n; i++){
          V(i,k)=w(i,0);
         } 
   //Solve upper hessenberg least squares problem via Givens rotations
   //Compute previous Givens rotations
    for (int i=0; i<k-1; i++){
       h(i,k-1)=c(i,0)*h(i,k-1)+s(i,0)*h(i+1,k-1);
       h(i+1,k-1)=-s(i,0)*h(i,k-1)+c(i,0)*h(i+1,k-1);
     }  
     //Compute next Givens rotations
     c.reshape(k,1);
     s.reshape(k,1); 
     bb.reshape(k+1,1);
     double l = sqrt(h(k-1,k-1)*h(k-1,k-1)+h(k,k-1)*h(k,k-1));
     c(k-1,0)=h(k-1,k-1)/l;
     s(k-1,0)=h(k,k-1)/l;
     std::cout <<" h(k,k-1)= " << h(k,k-1) << std::endl;
     // Givens rotation on h and bb
     h(k-1,k-1)=l;
     h(k,k-1)=0;
     bb(k-1,0)=c(k-1,0)*bb(k-1,0);
     bb(k,0)=-s(k-1,0)*bb(k-1,0);

    //Determine residual    
    resid = fabs(bb(k,0));

    std::cout << "resid = " << resid << std::endl;
    k=k+1;
  } 
   //Extract upper triangular square matrix
   bb.reshape(h.numRows()-1 ,1);

   //Solve linear system
   int info;

   Teuchos::LAPACK<int, double> lapack;
   lapack.TRTRS('U', 'N', 'N', h.numRows()-1, 1, h.values(), h.stride(), bb.values(), bb.stride(),&info); 
   //Found y=Mx
   for (int i=0; i<k-1; i++){
      bb(i,0)=bb(i,0)/D(i,0);
   }

   V.reshape(n,k-1);
   Teuchos::SerialDenseMatrix<int, double> ans(X);
   ans.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, V, bb, 1.0);
   std::cout << "ans= " << ans << std::endl;

   std::cout << "h= " << h << std::endl;



return 0;
}