00001
00037 #ifndef _APRX_H_
00038 #define _APRX_H_
00039
00040 #include <mc/mc_base.h>
00041
00042 #include <mc/gem.h>
00043 #include <mc/bam.h>
00044 #include <mc/whb.h>
00045
00046 #include <mc/re.h>
00047 #include <mc/node.h>
00048 #include <mc/bnode.h>
00049
00056 struct sAprx {
00057
00059 Vmem *vmem;
00061 int iMadeVmem;
00062
00064 Gem *gm;
00066 PDE *pde;
00067
00069 Vset *lnkG;
00070
00072 Re *re[MAXV];
00074 Re *reB[MAXV];
00075
00077 int numV;
00079 int numE;
00081 int numF;
00083 int numS;
00084
00086 int *I_V2E;
00088 int *I_V2F;
00090 int *I_V2S;
00092 int *I_E2F;
00094 int *I_E2S;
00096 int *I_F2S;
00097
00099 Bnode *B;
00100
00102 int numBV[MAXV];
00104 int *ibv[MAXV];
00106 double *bv[MAXV];
00107
00109 Bvec *wev;
00111 double gerror;
00112
00114 int numB;
00116 Bmat *P;
00117
00119 int order;
00120
00121 };
00122
00130 typedef struct sAprx Aprx;
00131
00132
00133
00134
00135
00136
00137
00138 #if !defined(VINLINE_APRX)
00139 #else
00140 #endif
00141
00152 Aprx* Aprx_ctor(Vmem *vmem, Gem *tgm, PDE *tpde);
00153
00165 Aprx* Aprx_ctor2(Vmem *vmem, Gem *tgm, PDE *tpde,int order);
00166
00174 void Aprx_dtor(Aprx **thee);
00175
00184 int Aprx_numB(Aprx *thee);
00185
00194 Bmat *Aprx_P(Aprx *thee);
00195
00217 int Aprx_read(Aprx *thee, int key, Vio *sock);
00218
00227 void Aprx_reset(Aprx *thee);
00228
00237 void Aprx_initNodes(Aprx *thee);
00238
00248 void Aprx_buildNodes(Aprx *thee);
00249
00261 void Aprx_evalTrace(Aprx *thee, Bvec *Wud, Bvec *Wui, Bvec *Wut);
00262
00298 int Aprx_refine(Aprx *thee, int rkey, int bkey, int pkey);
00299
00310 int Aprx_unRefine(Aprx *thee, int rkey, int pkey);
00311
00321 int Aprx_deform(Aprx *thee, Bvec *def);
00322
00335 void Aprx_buildProlong(Aprx *thee, int numRc, int numR, int pkey);
00336
00347 void Aprx_buildProlongBmat(Aprx *thee, Mat *p);
00348
00362 int Aprx_nodeCount(Aprx *thee, Re *re);
00363
00398 void Aprx_interact(Aprx *thee, Re *re, Re *reT,
00399 int *numR, int *numO, int *numOYR, int **IJA, int **IJAYR);
00400
00418 void Aprx_interactBlock(Aprx *thee, Re *re[MAXV], Re *reT[MAXV],
00419 MATsym sym[MAXV][MAXV], MATmirror mirror[MAXV][MAXV], MATformat frmt[MAXV][MAXV],
00420 int numR[MAXV], int numO[MAXV][MAXV], int *IJA[MAXV][MAXV]);
00421
00434 int Aprx_nodeIndex (Aprx *thee, Re *re, int i, int dim, int q);
00435
00449 int Aprx_nodeIndexTT (Aprx *thee, Re *re, TT *t, int idf);
00450
00463 void Aprx_buildBRC(Aprx *thee, Bmat *A, Bmat *M);
00464
00480 void Aprx_initSpace(Aprx *thee, int i, int r,
00481 double phi[], double phix[][3], double V[], double dV[][3]);
00482
00483
00484 #if 0
00485
00510 void Aprx_buildBasis(Aprx *thee, Re *re,
00511 int evalKey, int qp, int face, TT *t,
00512 double xq[], double phi[], double phix[][3],
00513 double U[], double dU[][3],
00514 Bvec *Wu, Bvec *Wud);
00515 #endif
00516
00530 void Aprx_initEmat(Aprx *thee, int bumpKey, TT *t, Emat *em);
00531
00580 void Aprx_quadEmat(Aprx *thee,
00581 int evalKey, int energyKey, int residKey, int tangKey, int massKey,
00582 int bumpKey,
00583 int face, TT *t, Emat *em,
00584 Bvec *Wu, Bvec *Wud);
00585
00624 void Aprx_fanEmat(Aprx *thee,
00625 int evalKey, int energyKey, int residKey, int tangKey, int massKey,
00626 int bumpKey,
00627 Emat *em,
00628 Bmat *A, Bmat *M, Bvec *F);
00629
00667 void Aprx_assemEmat(Aprx *thee, SS *sm,
00668 int evalKey, int energyKey, int residKey, int tangKey, int massKey,
00669 int bumpKey,
00670 Emat *em,
00671 Bvec *Wu, Bvec *Wud);
00672
00717 double Aprx_assem(Aprx *thee,
00718 int evalKey, int energyKey, int residKey, int tangKey, int massKey,
00719 int bumpKey,
00720 int ip[], double rp[],
00721 Bmat *A, Bmat *M, Bvec *F, Bvec *Wu, Bvec *Wud);
00722
00748 int Aprx_markRefine(Aprx *thee, int key, int color,
00749 int bkey, double elevel, Bvec *u, Bvec *ud, Bvec *f, Bvec *r);
00750
00785 int Aprx_estRefine(Aprx *thee, int key, int color,
00786 int bkey, double elevel, Bvec *u, Bvec *ud, Bvec *f, Bvec *r);
00787
00803 int Aprx_markRefineFixed (Aprx *thee, int num2ref, int color);
00804
00805
00820 void Aprx_evalFunc(Aprx *thee,
00821 Bvec *u, int block, int numPts, double *pts,
00822 double *vals, int *marks);
00823
00928 void Aprx_fe2fd(Aprx *thee, Bvec *u, int block,
00929 int nx, int ny, int nz,
00930 double x0, double y0, double z0, double dx, double dy, double dz,
00931 int derivs, int outTypeKey);
00932
00944 void Aprx_bndIntegral(Aprx *thee, Bvec *u, Bvec *ud, Bvec *ut);
00945
00958 void Aprx_admMass(Aprx *thee, Bvec *u, Bvec *ud, Bvec *ut, int block);
00959
00981 double Aprx_evalError(Aprx *thee, int pcolor, int key,
00982 Bvec *u, Bvec *ud, Bvec *ut);
00983
00997 double Aprx_evalH1(Aprx *thee, SS *sm, Bvec *u, Bvec *ud, Bvec *ut);
00998
01027 void Aprx_writeGEOM(Aprx *thee, Vio *sock,
01028 int defKey, int colKey, int chartKey, double gluVal, int fkey,
01029 Bvec *u, char *format);
01030
01042 void Aprx_writeSOL(Aprx *thee, Vio *sock, Bvec *u, char *format);
01043
01053 int Aprx_partSet(Aprx *thee, int pcolor);
01054
01063 int Aprx_partSmooth(Aprx *thee);
01064
01101 int Aprx_partInert(Aprx *thee, int pcolor,
01102 int numC, double *evec, simHelper *simH);
01103
01147 int Aprx_partSpect(Aprx *thee, int pcolor,
01148 int numC, double *evec, simHelper *simH, int *ford, int *rord,
01149 int general);
01150
01162 int Aprx_part(Aprx *thee, int pkey, int pwht, int ppow);
01163
01175 int Aprx_partOne(Aprx *thee, int pkey, int pwht, int pcolor, int poff);
01176
01177 #endif
01178
01179