• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

MOR_BasisOps.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 #ifndef MOR_BASISOPS_HPP
00007 #define MOR_BASISOPS_HPP
00008 
00009 #include "Epetra_MultiVector.h"
00010 
00011 class Epetra_LocalMap;
00012 class Epetra_Operator;
00013 
00014 namespace MOR {
00015 
00016 // Convenience functions for common reduced-order basis computations
00017 Epetra_LocalMap createComponentMap(const Epetra_MultiVector &projector);
00018 
00019 // result <- (basis * components) + (resultScaling * result)
00020 inline
00021 int expandAdd(const Epetra_MultiVector &basis, const Epetra_MultiVector &components, double resultScaling, Epetra_MultiVector &result)
00022 {
00023   return result.Multiply('N', 'N', 1.0, basis, components, resultScaling);
00024 }
00025 
00026 // result <- (basis * components) + result
00027 inline
00028 int expandAdd(const Epetra_MultiVector &basis, const Epetra_MultiVector &components, Epetra_MultiVector &result)
00029 {
00030   return expandAdd(basis, components, 1.0, result);
00031 }
00032 
00033 // result <- basis * components
00034 inline
00035 int expand(const Epetra_MultiVector &basis, const Epetra_MultiVector &components, Epetra_MultiVector &result)
00036 {
00037   return expandAdd(basis, components, 0.0, result);
00038 }
00039 
00040 // result <- (basis^T * vectors) + (resultScaling * result)
00041 inline
00042 int reduceAdd(const Epetra_MultiVector &basis, const Epetra_MultiVector &vectors, double resultScaling, Epetra_MultiVector &result)
00043 {
00044   return result.Multiply('T', 'N', 1.0, basis, vectors, resultScaling);
00045 }
00046 
00047 // result <- (basis^T * vectors) + result
00048 inline
00049 int reduceAdd(const Epetra_MultiVector &basis, const Epetra_MultiVector &vectors, Epetra_MultiVector &result)
00050 {
00051   return reduceAdd(basis, vectors, 1.0, result);
00052 }
00053 
00054 // result <- basis^T * vectors
00055 inline
00056 int reduce(const Epetra_MultiVector &basis, const Epetra_MultiVector &vectors, Epetra_MultiVector &result)
00057 {
00058   return reduceAdd(basis, vectors, 0.0, result);
00059 }
00060 
00061 // dual <- dual * (primal^T * dual)^{-1}
00062 void dualize(const Epetra_MultiVector &primal, Epetra_MultiVector &dual);
00063 
00064 // dual <- (metric * primal) * (primal^T * metric * primal)^{-1}
00065 void dualize(const Epetra_MultiVector &primal, const Epetra_Operator &metric, Epetra_MultiVector &result);
00066 
00067 } // namespace MOR
00068 
00069 #endif /* MOR_BASISOPS_HPP */

Generated on Wed Mar 26 2014 18:36:40 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1