SundanceFunctionSupportResolver.hpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 // ************************************************************************
00003 // 
00004 //                             Sundance
00005 //                 Copyright 2011 Sandia Corporation
00006 // 
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Kevin Long (kevin.long@ttu.edu)
00038 // 
00039 
00040 /* @HEADER@ */
00041 
00042 #ifndef SUNDANCE_FUNCTION_SUPPORT_RESOLVER_H
00043 #define SUNDANCE_FUNCTION_SUPPORT_RESOLVER_H
00044 
00045 #include "SundanceDefs.hpp"
00046 #include "SundanceExpr.hpp"
00047 #include "SundanceSet.hpp"
00048 #include "SundanceMap.hpp"
00049 #include "SundanceOrderedHandle.hpp"
00050 #include "SundanceCellFilterStub.hpp"
00051 
00052 
00053 namespace Sundance
00054 {
00055 
00056 class SumOfIntegrals;
00057 class SumOfBCs;
00058 class CommonFuncDataStub;
00059 
00060 /** */
00061 class FunctionSupportResolver
00062 {
00063 public:
00064   /** */
00065   FunctionSupportResolver(
00066     const Expr& eqns,
00067     const Expr& bcs,
00068     const Array<Expr>& vars,
00069     const Array<Expr>& unks,
00070     const Expr& unkParams,
00071     const Expr& params,
00072     const Array<Expr>& fixedFields,
00073     bool isVariational);
00074   
00075   
00076   /** \name Getting information about functions */
00077   //@{
00078   /** Returns the number of variational function blocks */
00079   int numVarBlocks() const {return varFuncs_.size();}
00080 
00081   /** Returns the number of unknown function blocks */
00082   int numUnkBlocks() const {return unkFuncs_.size();}
00083 
00084   /** Returns the number of unknown parameters */
00085   int numUnkParams() const {return unkParams_.size();}
00086 
00087   /** Returns the number of fixed parameters */
00088   int numFixedParams() const {return fixedParams_.size();}
00089 
00090   /** Returns the number of variational functions in this block */
00091   int numVars(int block) const {return varFuncs_[block].size();}
00092 
00093   /** Returns the number of unk functions in this block */
00094   int numUnks(int block) const {return unkFuncs_[block].size();}
00095 
00096   /** Returns the number of variational function IDs in this block.
00097    * This will differ from the number of variational functions in cases
00098    * where a vector field uses a single vector-valued basis rather
00099    * than scalar bases for each component. */
00100   int numVarIDs(int block) const {return varFuncData_[block].size();}
00101 
00102   /** Returns the number of unk function IDs in this block.
00103    * This will differ from the number of unknown functions in cases
00104    * where a vector field uses a single vector-valued basis rather
00105    * than scalar bases for each component. */ 
00106   int numUnkIDs(int block) const {return unkFuncData_[block].size();}
00107 
00108   /** Returns the data for the i-th variational function in block b */
00109   RCP<const CommonFuncDataStub> varFuncData(int b, int i) const {return varFuncData_[b][i];}
00110 
00111   /** Returns the data for the i-th unknown function in block b */
00112   RCP<const CommonFuncDataStub> unkFuncData(int b, int i) const {return unkFuncData_[b][i];}
00113 
00114   /** Returns the i-th unknown parameter */
00115   const Expr& unkParam(int i) const {return unkParams_[i];}
00116 
00117   /** Determine whether a given func ID is listed as a 
00118    * variational function in this equation set */
00119   bool hasVarID(int fid) const 
00120     {return varIDToBlockMap_.containsKey(fid);}
00121 
00122   /** Determine whether a given func ID is listed as a unk function 
00123    * in this equation set */
00124   bool hasUnkID(int fid) const 
00125     {return unkIDToBlockMap_.containsKey(fid);}
00126 
00127   /** Determine whether a given func ID is listed as a unk parameter 
00128    * in this equation set */
00129   bool hasUnkParamID(int fid) const 
00130     {return unkParamIDToReducedUnkParamIDMap_.containsKey(fid);}
00131 
00132   /** Determine whether a given func ID is listed as a fixed parameter 
00133    * in this equation set */
00134   bool hasFixedParamID(int fid) const 
00135     {return fixedParamIDToReducedFixedParamIDMap_.containsKey(fid);}
00136 
00137   /** get the block number for the variational function having the
00138    * specified unreduced funcID */
00139   int blockForVarID(int varID) const ;
00140 
00141   /** get the block number for the unknown function having the
00142    * specified unreduced funcID */
00143   int blockForUnkID(int unkID) const ;
00144   //@}
00145 
00146 
00147   /** \name Finding the functions that appear on regions */
00148   //@{
00149   /** Returns the variational functions that appear explicitly
00150    * on the d-th region */
00151   const Set<int>& varsOnRegion(int d) const 
00152     {return varsOnRegions_.get(regions_[d]);}
00153 
00154   /** Returns the unknown functions that appear explicitly on the
00155    * d-th region. */
00156   const Set<int>& unksOnRegion(int d) const 
00157     {return unksOnRegions_.get(regions_[d]);}
00158 
00159   /** Returns the variational functions that 
00160    * appear in BCs on the d-th region.
00161    * We can use this information to tag certain rows as BC rows */
00162   const Set<int>& bcVarsOnRegion(int d) const 
00163     {return bcVarsOnRegions_.get(regions_[d]);}
00164 
00165   /** Returns the unknown functions that appear in BCs on the d-th region.
00166    * We can use this information to tag certain columns as BC
00167    * columns in the event we're doing symmetrized BC application */
00168   const Set<int>& bcUnksOnRegion(int d) const 
00169     {return bcUnksOnRegions_.get(regions_[d]);}
00170 
00171   /** Returns the reduced variational functions that appear explicitly
00172    * on the d-th region */
00173   const Array<Set<int> >& reducedVarsOnRegion(const OrderedHandle<CellFilterStub>& r) const 
00174     {return reducedVarsOnRegions_[indexForRegion(r)];}
00175 
00176   /** Returns the reduced unknown functions that appear explicitly on the
00177    * d-th region. */
00178   const Array<Set<int> >& reducedUnksOnRegion(const OrderedHandle<CellFilterStub>& r) const 
00179     {return reducedUnksOnRegions_[indexForRegion(r)];}
00180   //@}
00181 
00182       
00183 
00184 
00185   /** \name Transforming between unreduced and reduced function IDs */
00186   //@{
00187   /** get the reduced ID for the variational function having the
00188    * specified unreduced funcID */
00189   int reducedVarID(int varID) const ;
00190 
00191   /** get the reduced ID for the unknown 
00192    * function having the given funcID */
00193   int reducedUnkID(int unkID) const ;
00194 
00195   /** get the reduced ID for the unk parameter
00196    * having the given funcID */
00197   int reducedUnkParamID(int unkID) const ;
00198 
00199   /** get the reduced ID for the fixed parameter
00200    * having the given funcID */
00201   int reducedFixedParamID(int unkID) const ;
00202 
00203   /** get the unreduced funcID for a variational function
00204    * as specified by a reduced ID and block index */
00205   int unreducedVarID(int block, int reducedVarID) const 
00206     {return unreducedVarID_[block][reducedVarID];}
00207 
00208   /** get the unreduced funcID for an unknown function
00209    * as specified by a reduced ID and block index */
00210   int unreducedUnkID(int block, int reducedUnkID) const 
00211     {return unreducedUnkID_[block][reducedUnkID];}
00212 
00213   /** get the unreduced funcID for an unknown parameter
00214    * as specified by a reduced ID  */
00215   int unreducedUnkParamID(int reducedUnkParamID) const 
00216     {return unreducedUnkParamID_[reducedUnkParamID];}
00217 
00218   /** get the unreduced funcID for a fixed parameter
00219    * as specified by a reduced ID */
00220   int unreducedFixedParamID(int reducedFixedParamID) const 
00221     {return unreducedFixedParamID_[reducedFixedParamID];}
00222 
00223   /** Return the map from fixed param ID to reduced fixed param ID */
00224   const Map<int, int>& fixedParamIDToReducedFixedParamIDMap() const
00225     {return fixedParamIDToReducedFixedParamIDMap_;}
00226   
00227   //@}
00228 
00229 
00230   /** \name Finding integration regions for the equation set */
00231   //@{
00232   /** Returns the number of regions on which pieces of the equation
00233    * or BCs are defined. */
00234   int numRegions() const {return regions_.size();}
00235       
00236   /** Returns the d-th region for this equation set */
00237   const RCP<CellFilterStub>& region(int d) const 
00238     {return regions_[d].ptr();}
00239 
00240   /** Returns the index of the given region */
00241   int indexForRegion(const OrderedHandle<CellFilterStub>& region) const ;
00242 
00243   /** Whether a region has BCs */
00244   bool isBCRegion(int d) const 
00245     {return bcVarsOnRegions_.containsKey(regions_[d]);}
00246 
00247 
00248   /** Return the set of regions on which the specified 
00249    * test func appears. */
00250   const Set<OrderedHandle<CellFilterStub> >& 
00251   regionsForTestFunc(int unreducedTestID) const ;
00252       
00253   /** Return the set of regions on which the specified 
00254    * unknown func appears */
00255   const Set<OrderedHandle<CellFilterStub> >& 
00256   regionsForUnkFunc(int unreducedUnkID) const ;
00257   //@}
00258 
00259 
00260   /**
00261    * Flatten a spectral expression into a list of its coefficients
00262    */
00263   Expr flattenSpectral(const Expr& input) const ;
00264   /**
00265    * Flatten a spectral expression into a list of its coefficients
00266    */
00267   Array<Expr> flattenSpectral(const Array<Expr>& input) const ;
00268 
00269   /** Whether essential BCs are present */
00270   bool hasBCs() const ;
00271 
00272   /** Access to integrals */
00273   const SumOfIntegrals* integralSum() const {return integralSum_;}
00274 
00275   /** Access to BCs */
00276   const SumOfBCs* bcSum() const {return bcSum_;}
00277 
00278   /** */
00279   const Array<Expr>& unks() const {return unkFuncs_;}
00280 
00281   /** */
00282   const Array<Expr>& vars() const {return varFuncs_;}
00283 
00284   /** */
00285   const Array<Expr>& fixedFields() const {return fixedFields_;}
00286 
00287   /** */
00288   const Expr& fixedParams() const {return fixedParams_;}
00289 
00290   /** */
00291   const Expr& unkParams() const {return unkParams_;}
00292 
00293   /** */
00294   const Set<int>& varFuncSet() const {return varFuncSet_;}
00295   /** */
00296   const Set<int>& unkFuncSet() const {return unkFuncSet_;}
00297   /** */
00298   const Set<int>& unkParamSet() const {return unkParamSet_;}
00299   /** */
00300   const Set<int>& fixedParamSet() const {return fixedParamSet_;}
00301   
00302       
00303 private:
00304 
00305   /** */
00306   Expr eqns_;
00307 
00308   /** */
00309   Expr bcs_;
00310 
00311   /** */
00312   const SumOfIntegrals* integralSum_;
00313 
00314   /** */
00315   const SumOfBCs* bcSum_;
00316 
00317   /** */
00318   Set<int> varFuncSet_;
00319 
00320   /** */
00321   Set<int> unkFuncSet_;
00322 
00323   /** */
00324   Set<int> unkParamSet_;
00325 
00326   /** */
00327   Set<int> fixedParamSet_;
00328 
00329   /** */
00330   Array<OrderedHandle<CellFilterStub> > regions_;
00331 
00332   /** */
00333   Map<OrderedHandle<CellFilterStub>, int> regionToIndexMap_;
00334 
00335   /** */
00336   Map<OrderedHandle<CellFilterStub>, Set<int> > varsOnRegions_;
00337 
00338   /** */
00339   Map<OrderedHandle<CellFilterStub>, Set<int> > unksOnRegions_;
00340 
00341   /** */
00342   Array<Array<Set<int> > > reducedVarsOnRegions_;
00343 
00344   /** */
00345   Array<Array<Set<int> > > reducedUnksOnRegions_;
00346 
00347   /** */
00348   Map<OrderedHandle<CellFilterStub>, Set<int> > bcVarsOnRegions_;
00349 
00350   /** */
00351   Map<OrderedHandle<CellFilterStub>, Set<int> > bcUnksOnRegions_;
00352 
00353   /** */
00354   Map<int, Set<OrderedHandle<CellFilterStub> > > testToRegionsMap_;
00355 
00356   /** */
00357   Map<int, Set<OrderedHandle<CellFilterStub> > > unkToRegionsMap_;
00358 
00359   /** var function data for this equation set */
00360   Array<Array<RCP<const CommonFuncDataStub> > > varFuncData_;
00361 
00362   /** unknown function data for this equation set */
00363   Array<Array<RCP<const CommonFuncDataStub> > > unkFuncData_;
00364 
00365   /** var functions for this equation set */
00366   Array<Expr> varFuncs_;
00367 
00368   /** unknown functions for this equation set */
00369   Array<Expr> unkFuncs_;
00370 
00371   /** fixed functions for this equation set */
00372   Array<Expr> fixedFields_;
00373 
00374   /** The point in function space about which the equations
00375    * are linearized */
00376   Array<Expr> unkLinearizationPts_;
00377 
00378   /** unknown parameters for this equation set */
00379   Expr unkParams_;
00380 
00381   /** fixed parameters for this equation set */
00382   Expr fixedParams_;
00383 
00384   /** map from variational function funcID to that function's
00385    * position in list of var functions */
00386   Array<Map<int, int> > varIDToReducedIDMap_;
00387 
00388   /** map from unknown function funcID to that function's
00389    * position in list of unk functions */
00390   Array<Map<int, int> > unkIDToReducedIDMap_;
00391 
00392   /** map from unknown param funcID to that param's
00393    * position in list of unk params */
00394   Map<int, int> unkParamIDToReducedUnkParamIDMap_;
00395 
00396   /** map from fixed param funcID to that param's
00397    * position in list of fixed params */
00398   Map<int, int> fixedParamIDToReducedFixedParamIDMap_;
00399 
00400   /** map from variational function funcID to that function's
00401    * position in list of var blocks */
00402   Map<int, int> varIDToBlockMap_;
00403 
00404   /** map from unknown function funcID to that function's
00405    * position in list of unk blocks */
00406   Map<int, int> unkIDToBlockMap_;
00407 
00408   /** Map from (block, unreduced var ID) to reduced ID */
00409   Array<Array<int> > reducedVarID_;
00410 
00411   /** Map from (block, unreduced unk ID) to reduced ID */
00412   Array<Array<int> > reducedUnkID_;
00413 
00414   /** Map from unreduced unk ID to reduced ID */
00415   Array<int> reducedUnkParamID_;
00416 
00417   /** Map from unreduced fixed ID to reduced ID */
00418   Array<int> reducedFixedParamID_;
00419 
00420   /** Map from (block, reduced varID) to unreduced varID */
00421   Array<Array<int> > unreducedVarID_;
00422 
00423   /** Map from (block, reduced unkID) to unreduced unkID */
00424   Array<Array<int> > unreducedUnkID_;
00425 
00426   /** Map from reduced unkParamID to unreduced unkParamID */
00427   Array<int> unreducedUnkParamID_;
00428 
00429   /** Map from reduced fixedParamID to unreduced fixedParamID */
00430   Array<int> unreducedFixedParamID_;
00431 
00432   /** Flag indicating whether this equation set is nonlinear */
00433   bool isNonlinear_;
00434       
00435   /** Flag indicating whether this equation set is 
00436    * a variational problem */
00437   bool isVariationalProblem_;
00438 };
00439 
00440 /** */
00441 RCP<const CommonFuncDataStub> getSharedFunctionData(const FuncElementBase* f);
00442 
00443 }
00444  
00445 #endif

Site Contact