00001 #ifndef WORKSET_H
00002 #define WORKSET_H
00003
00004 #include "Cell.h"
00005 #include "CellType.h"
00006 #include "DenseSerialVector.h"
00007 #include "TSFArray.h"
00008 #include "QuadraturePoints.h"
00009 #include "Expr.h"
00010 #include "QuadBasisSpecifier.h"
00011 #include "QuadExprSpecifier.h"
00012 #include "TSFHashtable.h"
00013 #include "MultiIndex.h"
00014 #include "DistributedDOFMap.h"
00015 #include "QuadratureFamily.h"
00016 #include "LocalMatrix.h"
00017 #include "WeakForm.h"
00018 #include "UserDefinedFunction.h"
00019 #include "TSFTimeMonitor.h"
00020
00021 namespace Sundance
00022 {
00023 using TSF::DenseSerialVector;
00024
00025
00026
00027
00028 class WorkSet
00029 {
00030 public:
00031 WorkSet(int capacity, const Cell& prototype);
00032
00033
00034 void setSize(int nCells) {nCells_ = nCells;}
00035
00036
00037 int nCells() const {return nCells_;}
00038
00039
00040 int capacity() const {return capacity_;}
00041
00042
00043 int cellDim() const {return cellDim_;}
00044
00045
00046 int spatialDim() const {return spatialDim_;}
00047
00048
00049 int cellIndex(int c) const {return cellIndex_[c];}
00050
00051
00052 const TSFArray<Cell> cells() const {return cells_;}
00053
00054
00055
00056
00057
00058
00059 const CellTopologyCode& integrationCellType(int derivOrder) const ;
00060
00061
00062 const CellTopologyCode& cellType() const {return cellType_;}
00063
00064
00065 const Cell& getCell(int c) const {return cells_[c];}
00066
00067
00068 void addCell(const Cell& cell, int& cellCount);
00069
00070
00071 void resetCaches();
00072
00073
00074 void buildLocalMatrix(const WeakForm& weakForm,
00075 TSFArray<LocalMatrix>& matrix) ;
00076
00077 void buildLocalVector(const WeakForm& wekaForm,
00078 TSFArray<DenseSerialVector>& vector) ;
00079
00080
00081 double integrate(const Expr& integrand, const QuadratureFamily& quad) ;
00082
00083
00084
00085 TSFSmartPtr<DenseSerialVector>
00086 getPhysBasisValues(const BasisFamily& basis,
00087 const QuadratureFamily& quad,
00088 int derivOrder) const ;
00089
00090 TSFSmartPtr<DenseSerialVector>
00091 getExprValues(const Expr& expr,
00092 const QuadratureFamily& quad) const ;
00093
00094
00095 void getCoordExprValues(const QuadratureFamily& qp,
00096 int direction,
00097 DenseSerialVector& result) const ;
00098
00099
00100
00101 void getUserFuncExprValues(const QuadratureFamily& qp,
00102 const TSFSmartPtr<UserDefinedFunction>& func,
00103 DenseSerialVector& result) const ;
00104
00105
00106
00107 void getCellDiameterExprValues(const QuadratureFamily& qp,
00108 DenseSerialVector& result) const ;
00109
00110
00111 void getCellNormalExprValues(const QuadratureFamily& qp,
00112 int direction,
00113 DenseSerialVector& result) const ;
00114
00115
00116 void evaluateDiscreteFunc(const TSFSmartPtr<DOFMapBase>& map,
00117 const CellSet& domain,
00118 const TSFVector& vector,
00119 const BasisFamily& basis,
00120 int funcIndex,
00121 const MultiIndex& deriv,
00122 const QuadratureFamily& quadFamily,
00123 DenseSerialVector& result) const ;
00124
00125
00126 TSFSmartPtr<DenseSerialVector> getQuadWeights(const QuadratureFamily& q) const ;
00127
00128 int getQuadSize(const QuadratureFamily& q) const ;
00129
00130 void setupQuadrature(const QuadratureFamily& q) const ;
00131
00132 inline const int* cellIndex() const {return &(cellIndex_[0]);}
00133
00134 inline const double* detJ() const {return &(detJ_[0]);}
00135
00136 inline const double& detJ(int c, int numQuadPoints = 1, int quadPoint=0) const {return detJ_[c*numQuadPoints+quadPoint];}
00137
00138 static void verboseOn() {verbose_ = true;}
00139 static void verboseOff() {verbose_ = false;}
00140
00141 static bool& paranoidChecking() {static bool rtn = false; return rtn;}
00142 private:
00143
00144
00145 TSFSmartPtr<DenseSerialVector>
00146 getRefBasisValues(const BasisFamily& basis,
00147 const QuadratureFamily& quad,
00148 int derivOrder) const ;
00149
00150 void addMaximalCell(const Cell& cell, int cellCount);
00151
00152 void addEmbeddedCell(const Cell& cell, int cellCount);
00153
00154 void setJacobianValues(const QuadratureFamily& quadFamily);
00155
00156 void getJacobian(const Cell& cell, int cellCount, const QuadraturePoints& quadPoints);
00157
00158 void getConstantJacobian(const Cell& cell, int cellCount);
00159
00160
00161 TSFSmartPtr<DenseSerialVector> getPhysicalQuadPoints(const QuadratureFamily& quad) const ;
00162
00163
00164 TSFSmartPtr<DenseSerialVector> getReferenceQuadPoints(const QuadratureFamily& quad) const ;
00165
00166
00167 TSFSmartPtr<DenseSerialVector> getParentQuadPoints(const QuadratureFamily& quad) const ;
00168
00169
00170 void computeCoords(const QuadratureFamily& quadFamily) const ;
00171
00172
00173 void computeCellDiameters(const QuadratureFamily& quadFamily) const ;
00174
00175
00176 void computeCellNormals(const QuadratureFamily& quadFamily) const ;
00177
00178 inline double* J() {return &(J_[0]);}
00179
00180 inline const double* invJ() const {return &(invJ_[0]);}
00181
00182 inline double* invJ() {return &(invJ_[0]);}
00183
00184
00185
00186 int capacity_;
00187
00188
00189 int nCells_;
00190
00191
00192 int spatialDim_;
00193
00194
00195 int cellDim_;
00196
00197
00198 CellTopologyCode cellType_;
00199
00200
00201 CellTopologyCode maximalCellType_;
00202
00203
00204 int nodesPerCell_;
00205
00206
00207
00208 bool isAffine_;
00209
00210
00211 TSFArray<int> cellIndex_;
00212
00213
00214 TSFArray<Cell> cells_;
00215
00216
00217
00218 DenseSerialVector detJ_;
00219
00220
00221
00222 DenseSerialVector J_;
00223
00224
00225
00226 DenseSerialVector invJ_;
00227
00228
00229 DenseSerialVector jTmp_;
00230
00231 TSFArray<int> iPiv_;
00232
00233
00234 mutable TSFHashtable<QuadratureFamily,
00235 TSFSmartPtr<DenseSerialVector> > refQuadPts_;
00236
00237
00238 mutable TSFHashtable<QuadratureFamily,
00239 TSFSmartPtr<DenseSerialVector> > physQuadPts_;
00240
00241
00242 mutable TSFHashtable<QuadratureFamily,
00243 TSFSmartPtr<DenseSerialVector> > parentQuadPts_;
00244
00245
00246 mutable TSFHashtable<QuadratureFamily, int> nQuadPoints_;
00247
00248
00249 mutable TSFHashtable<QuadratureFamily,
00250 TSFSmartPtr<DenseSerialVector> > quadWeights_;
00251
00252
00253
00254 mutable TSFHashtable<QuadBasisSpecifier,
00255 TSFSmartPtr<DenseSerialVector> > refBasisValues_;
00256
00257
00258
00259 mutable TSFHashtable<QuadBasisSpecifier,
00260 TSFSmartPtr<DenseSerialVector> > physBasisValues_;
00261
00262
00263 mutable TSFHashtable<QuadExprSpecifier,
00264 TSFSmartPtr<DenseSerialVector> > exprValues_;
00265
00266
00267 mutable TSFHashtable<QuadratureFamily,
00268 TSFArray<TSFSmartPtr<DenseSerialVector> > > coordValues_;
00269
00270
00271 mutable TSFHashtable<QuadratureFamily,
00272 TSFArray<TSFSmartPtr<DenseSerialVector> > > cellNormalValues_;
00273
00274
00275 mutable TSFHashtable<QuadratureFamily,
00276 TSFSmartPtr<DenseSerialVector> > cellDiameterValues_;
00277
00278
00279
00280 mutable bool needsPhysQuadPts_;
00281
00282
00283
00284 mutable bool needsPhysBasisValues_;
00285
00286
00287
00288 mutable bool needsCoordValues_;
00289
00290
00291
00292
00293 mutable bool needsCellNormalValues_;
00294
00295
00296
00297
00298 mutable bool needsCellDiameterValues_;
00299
00300
00301 mutable bool needsJacobianValues_;
00302
00303
00304 mutable bool needsParentQuadPts_;
00305
00306
00307 void printRefQuadPoints(const QuadratureFamily& quad, int c) const ;
00308
00309
00310 void printPhysQuadPoints(const QuadratureFamily& quad, int c) const ;
00311
00312
00313 static TSFTimer& localMatrixTimer();
00314
00315
00316 static TSFTimer& localVectorTimer();
00317
00318
00319 static TSFTimer& integrationTimer();
00320
00321
00322 static TSFTimer& workSetAssemblyTimer();
00323
00324
00325 static TSFTimer& basisTransformationTimer();
00326
00327
00328 static TSFTimer& discreteFuncEvalTimer();
00329
00330
00331 static bool verbose_;
00332 };
00333 }
00334
00335 namespace TSF
00336 {
00337 inline string toString(const Sundance::WorkSet& ws)
00338 {
00339 return "WorkSet(" + ws.cellType().toString() + ")";
00340 }
00341 }
00342
00343 #endif