00001 #ifndef EXPR_H
00002 #define EXPR_H
00003
00004 #include "SundanceDefs.h"
00005
00006 #include "TSFArray.h"
00007 #include "BasisFamily.h"
00008 #include "DenseSerialVector.h"
00009 #include "TSFVector.h"
00010 #include "ExprHandle.h"
00011 #include "ExprArray.h"
00012 #include "LogicalExpr.h"
00013 #include "ExprBase.h"
00014 #include "Defaults.h"
00015
00016 namespace Sundance
00017 {
00018
00019 using namespace TSF;
00020 using std::string;
00021
00022 using std::ostream;
00023
00024
00025 class MultiIndex;
00026 class QuadratureFamily;
00027 class WorkSet;
00028 class Cell;
00029 class CellSet;
00030 class AbstractFunctionSpace;
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 class Expr
00044 {
00045 public:
00046
00047
00048
00049 Expr();
00050
00051 Expr(ExprBase* ptr);
00052
00053 Expr(const double& value);
00054
00055
00056
00057 ~Expr();
00058 Expr(const Expr& other);
00059 const Expr& operator=(const Expr& other);
00060
00061
00062
00063
00064 Expr operator-() const ;
00065
00066 const Expr& operator+=(const Expr& other) ;
00067
00068 const Expr& operator-=(const Expr& other) ;
00069
00070 const Expr& operator*=(const Expr& other) ;
00071
00072 const Expr& operator/=(const Expr& other) ;
00073
00074 Expr pow(const double& p) const ;
00075
00076
00077
00078
00079
00080 const Expr& operator[](int i) const ;
00081
00082 Expr& operator[](int i) ;
00083
00084 int length() const ;
00085
00086 int size() const ;
00087
00088 Expr append(const Expr& other);
00089
00090 Expr flatten() const ;
00091
00092
00093
00094
00095
00096 Expr variation(const Expr& u) const ;
00097
00098
00099
00100 Expr variation() const ;
00101
00102
00103 Expr linearization(const Expr& u, const Expr& u0) const ;
00104
00105
00106
00107 Expr directSensitivity(const Expr& u, const Expr& u0) const ;
00108
00109
00110
00111 Expr sensitivityEquation(const Expr& u, const Expr& alpha,
00112 const Expr& alpha0) const ;
00113
00114
00115
00116 Expr sensitivity(const Expr& alpha) const ;
00117
00118
00119
00120 Expr differential() const ;
00121
00122
00123 Expr substitute(const Expr& u, const Expr& u0) const ;
00124
00125
00126 Expr differential(const Expr& u, const Expr& du) const ;
00127
00128
00129 void setDifferential(const Expr& du) ;
00130
00131
00132 void setVariation(const Expr& v) ;
00133
00134
00135 void setFunctionValue(const Expr& u0) ;
00136
00137
00138
00139
00140
00141 LogicalExpr operator==(const Expr& other) const ;
00142
00143 LogicalExpr operator!=(const Expr& other) const ;
00144
00145 LogicalExpr operator<(const Expr& other) const ;
00146
00147 LogicalExpr operator>(const Expr& other) const ;
00148
00149 LogicalExpr operator<=(const Expr& other) const ;
00150
00151 LogicalExpr operator>=(const Expr& other) const ;
00152
00153
00154
00155
00156 ExprValue average(const Cell& cell) const ;
00157
00158 ExprValue average(const Mesh& mesh, const CellSet& cellSet) const ;
00159
00160 double probeAtMeshPoint(int localPointIndex) const ;
00161
00162
00163 void evaluate(const AbstractFunctionSpace& targetSpace, const TSFArray<Cell>& cells,
00164 const TSFArray<int>& cellIndices,
00165 const TSFArray<int>& dofIndices,
00166 const TSFArray<Point>& x,
00167 DenseSerialVector& values) const ;
00168
00169
00170 double evaluateFunctional(const Mesh& mesh,
00171 const Expr& u, const Expr& u0) const ;
00172
00173
00174 double evaluateFunctional(const Mesh& mesh) const ;
00175
00176
00177
00178
00179
00180 const Expr& getRegionalExpr(const string& region) const ;
00181
00182
00183 bool isDefinedOnRegion(const string& region) const ;
00184
00185
00186 void setRegionalExpr(const string& region, const Expr& expr);
00187
00188
00189
00190
00191
00192
00193 double integral(const Mesh& mesh,
00194 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00195
00196
00197 double integral(const Mesh& mesh,
00198 const CellSet& cellSet,
00199 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00200
00201
00202
00203 double integral(const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00204
00205
00206
00207 double integral(const CellSet& cellSet,
00208 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00209
00210
00211
00212 double norm(int p,
00213 const Mesh& mesh,
00214 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00215
00216
00217 double norm(int p,
00218 const Mesh& mesh,
00219 const CellSet& cellSet,
00220 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00221
00222
00223
00224 double norm(int p,
00225 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00226
00227
00228
00229 double norm(int p,
00230 const CellSet& cellSet,
00231 const QuadratureFamily& quad = Defaults::quadratureFamily()) const ;
00232
00233
00234
00235
00236
00237 double quickNorm() const ;
00238
00239
00240 double maxNorm() const ;
00241
00242
00243
00244
00245
00246 string toString(bool paren=false) const ;
00247
00248 string fullForm() const ;
00249
00250 XMLObject toXML() const ;
00251
00252 void matlabDump(ostream& os) const ;
00253
00254 void printDiscreteFunctions(ostream& os) const ;
00255
00256
00257
00258
00259
00260 void differential(const Expr& u, const Expr& du, Expr& result) const ;
00261
00262
00263 Expr derivative(const MultiIndex& d) const ;
00264
00265 const MultiIndex& multiIndex() const ;
00266
00267
00268 const Expr& assignUsingLHSPtr(const Expr& other);
00269
00270 const Expr& assignUsingRHSPtr(const Expr& other);
00271
00272
00273 double value() const ;
00274 void setParameterValue(const double& alpha);
00275
00276 const string& name() const ;
00277
00278
00279
00280
00281 Expr hold() const ;
00282
00283
00284 void evaluate(const WorkSet& workSet,
00285 const QuadratureFamily& quadFamily,
00286 DenseSerialVector& result) const ;
00287
00288
00289 void watchOn() ;
00290
00291 void watchOff() ;
00292
00293
00294 bool watching() const ;
00295
00296
00297 void getVector(TSFVector& vector) const ;
00298
00299
00300 void setVector(const TSFVector& vector) ;
00301
00302
00303 void readValues(const string& filename) ;
00304
00305
00306 void getDOFMap(TSFSmartPtr<DOFMapBase>& map) const ;
00307
00308
00309 int getReducedIndex() const ;
00310
00311
00312
00313 bool equals(const Expr& other) const ;
00314
00315 bool lessThan(const Expr& other) const ;
00316
00317 bool greaterThan(const Expr& other) const ;
00318
00319
00320 int sortPriority() const ;
00321
00322
00323 bool isNull() const {return handle()==0;}
00324
00325 bool isZero() const ;
00326
00327 bool isConstant() const ;
00328
00329 bool isSumExpr() const ;
00330
00331 bool isProductExpr() const ;
00332
00333 bool isTermListExpr() const ;
00334
00335 bool isFuncExpr() const ;
00336
00337 bool isDerivative() const ;
00338
00339 bool isDiffOp() const ;
00340
00341 bool isVariational() const ;
00342
00343 bool isUnknown() const ;
00344
00345 bool isCoordExpr() const ;
00346
00347 bool isDiscreteFunction() const ;
00348
00349 bool isUserFuncExpr() const ;
00350
00351 bool isParameterExpr() const ;
00352
00353 bool isHoldExpr() const ;
00354
00355 bool isRegionalExpr() const ;
00356
00357 bool isCellDiameterExpr() const ;
00358
00359 bool isCellNormalExpr() const ;
00360
00361 bool isIntegralExpr() const ;
00362
00363 bool isTestParameter() const ;
00364
00365 bool isUnknownParameter() const ;
00366
00367 bool isSpatiallyConstant() const ;
00368
00369
00370 bool containsUnknown() const ;
00371
00372 bool containsTest() const ;
00373
00374 bool containsUnknownParameter() const ;
00375
00376 bool containsTestParameter() const ;
00377
00378
00379 bool hasChildren() const ;
00380
00381
00382 void getChildren(ExprArray& children) const ;
00383
00384
00385
00386 int funcID() const ;
00387
00388
00389
00390 double constant() const ;
00391
00392
00393 const ExprBase* ptr() const {return handle()->ptr();}
00394
00395
00396
00397 Expr extractPrefactor(double& prefactor) const ;
00398
00399
00400 Expr clone() const ;
00401
00402
00403 bool lhsAppearsOnRHS(const Expr& lhs) const;
00404
00405 void replaceWithCloneOfAssignmentLHS(const Expr& lhs,
00406 const Expr& lhsClone);
00407
00408 void print(ostream& os, bool paren=false) const ;
00409
00410
00411
00412 BasisFamily getBasis() const ;
00413
00414 int countMonomials() const ;
00415
00416 ExprArray getMonomials() const ;
00417
00418 void getMonomials(ExprArray& monomials, int& offset) const ;
00419
00420
00421
00422
00423 bool getValidWeakForm(Expr& coeff, Expr& var, Expr& unk) const ;
00424
00425 bool getMesh(Mesh& mesh) const ;
00426
00427 void getDomain(TSFNonDupArray<CellSet>& domain) const ;
00428
00429
00430
00431
00432 friend Expr List(const Expr& a);
00433
00434 friend Expr List(const Expr& a, const Expr& b);
00435
00436 friend Expr List(const Expr& a, const Expr& b, const Expr& c);
00437
00438 friend Expr List(const Expr& a, const Expr& b, const Expr& c, const Expr& d);
00439
00440 int hashCode() const ;
00441
00442 static Expr& dummy();
00443
00444
00445 private:
00446 void considerSuicide();
00447 bool readyToDie() const ;
00448 void initializeRefCount();
00449 void initializeHandle(ExprHandle* handlePtr);
00450 void decrementRefCount() {(*refCount_)--;}
00451 void incrementRefCount() {(*refCount_)++;}
00452 ExprBase* ptr() {return handle()->ptr();}
00453 void nullCheck(const string& funcName) const ;
00454 ExprHandle* handle() {return *ptrToHandle_;}
00455 const ExprHandle* handle() const {return *ptrToHandle_;}
00456
00457
00458 ExprHandle** ptrToHandle_;
00459 int* refCount_;
00460
00461 };
00462
00463
00464 inline string toString(const Expr& expr)
00465 {
00466 return expr.toString();
00467 }
00468
00469
00470 inline ostream& operator<<(ostream& os, const Expr& expr)
00471 {
00472 expr.print(os);
00473 return os;
00474 }
00475
00476 inline int hashCode(const Expr& expr)
00477 {
00478 return expr.hashCode();
00479 }
00480
00481
00482
00483
00484
00485 Expr join(const Expr& e1, const Expr& e2);
00486
00487
00488
00489
00490
00491
00492 Expr operator+(const Expr& a, const Expr& b);
00493
00494
00495
00496
00497
00498 Expr operator-(const Expr& a, const Expr& b);
00499
00500
00501
00502
00503
00504 Expr operator*(const Expr& a, const Expr& b);
00505
00506
00507
00508
00509
00510 Expr operator/(const Expr& a, const Expr& b);
00511
00512
00513
00514
00515
00516 Expr pow(const Expr& expr, const double& p);
00517
00518
00519 Expr Scalar();
00520
00521
00522
00523
00524 Expr hold(const Expr& expr);
00525 }
00526 #endif