/* ****************************************
 * ( d[uvw]/d[xyz]*phi_j, phi_i )
 * return static float[9] managed by user.
 *  */
void
phgQuadGradDofBasDotBas(ELEMENT *e, DOF *gradu, DOF *u, int m, int n,
			int order, FLOAT *value)
{
    int i;
    const FLOAT *g1, *g2, *gu, *w, *lambda;
    FLOAT d[9], *pd, vol;
    QUAD *quad;

    assert(gradu->dim == 9 && u->dim == 3);

    if (order < 0)
	order = DofTypeOrder(gradu, e) +
	       BasisOrder(u, e, m) + BasisOrder(u, e, n) ;
    if (order < 0)
	order = 0;
    quad = phgQuadGetQuad3D(order);

    gu = phgQuadGetDofValues(e, gradu, quad);
    g1 = phgQuadGetBasisValues(e, u, m, quad);
    g2 = phgQuadGetBasisValues(e, u, n, quad);
    lambda = quad->points;
    w = quad->weights;
    for (i = 0; i < 9; i++)
	d[i] = 0;
    for (i = 0; i < quad->npoints; i++) {
	pd = d;
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* du/dx phi_m phi_n */
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* du/dy phi_m phi_n */
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* du/dz phi_m phi_n */

	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* dv/dx phi_m phi_n */
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* dv/dy phi_m phi_n */
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* dv/dz phi_m phi_n */

	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* dw/dx phi_m phi_n */
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* dw/dy phi_m phi_n */
	(*pd++) += (*(gu++)) * (*g1) * (*g2) * (*w);	/* dw/dz phi_m phi_n */

	g1++;
	g2++;
	w++;
	lambda += Dim + 1;
    }

    vol = phgGeomGetVolume(u->g, e);
    for (i = 0; i < 9; i++)
	value[i] = d[i] * vol;

    return;
}
/* *******************************************************
 * ( (u.\grad) phi_j, phi_i )
 * return 1 value
 *  */
FLOAT
phgQuadDofDotGradBasBas(ELEMENT *e, DOF *u, DOF *v, int m, int n, int order)
{
    int i;
    const FLOAT *g1, *g2, *gu, *w, *lambda;
    FLOAT d, d0;
    QUAD *quad;

    assert(u->dim == 3);

    if (order < 0)
	order = DofTypeOrder(u, e) +
	      BasisOrder(v, e, m) - 1 + BasisOrder(v, e, n);
    if (order < 0)
	order = 0;
    quad = phgQuadGetQuad3D(order);

    gu = phgQuadGetDofValues(e, u, quad);
    g1 = phgQuadGetBasisGradient(e, v, m, quad);
    g2 = phgQuadGetBasisValues(e, v, n, quad);
    d = 0.;
    lambda = quad->points;
    w = quad->weights;
    for (i = 0; i < quad->npoints; i++) {
	d0 = 0.;
	d0 += (*(gu++)) * (*(g1++)) * (*(g2));	/* u dphi_m/dx phi_n */
	d0 += (*(gu++)) * (*(g1++)) * (*(g2));	/* v dphi_m/dy phi_n */
	d0 += (*(gu++)) * (*(g1++)) * (*(g2));	/* w dphi_m/dz phi_n */
	g2++;
	d += d0 * (*(w++));
	lambda += Dim + 1;
    }

    return d * phgGeomGetVolume(u->g, e);
}
/* ***************************
 * ( (u.\grad) u, phi_i )
 * return 3 values.
 *  */
void
phgQuadDofDotGradDofBas(ELEMENT *e, DOF *u, DOF *gradu,
			int n, int order, FLOAT *values)
{
    int i;
    const FLOAT *g2, *gu, *ggu, *w, *lambda;
    FLOAT d1, d2, d3, dd1, dd2, dd3, vol;
    QUAD *quad;

    assert(gradu->dim == 9 && u->dim == 3);

    if (order < 0)
	order = DofTypeOrder(u, e) +
	       DofTypeOrder(gradu, e) + BasisOrder(u, e, n);
    if (order < 0)
	order = 0;
    quad = phgQuadGetQuad3D(order);

    gu = phgQuadGetDofValues(e, u, quad);
    ggu = phgQuadGetDofValues(e, gradu, quad);
    g2 = phgQuadGetBasisValues(e, u, n, quad);
    dd1 = dd2 = dd3 = 0.;
    lambda = quad->points;
    w = quad->weights;
    for (i = 0; i < quad->npoints; i++) {
	d1 = d2 = d3 = 0.;

	d1 += *(gu) * (*(ggu++)) * (*(g2));	/* u du/dx phi_n */
	d1 += *(gu + 1) * (*(ggu++)) * (*(g2));	/* v du/dy phi_n */
	d1 += *(gu + 2) * (*(ggu++)) * (*(g2));	/* w du/dz phi_n */

	d2 += *(gu) * (*(ggu++)) * (*(g2));	/* u dv/dx phi_n */
	d2 += *(gu + 1) * (*(ggu++)) * (*(g2));	/* v dv/dy phi_n */
	d2 += *(gu + 2) * (*(ggu++)) * (*(g2));	/* w dv/dz phi_n */

	d3 += *(gu) * (*(ggu++)) * (*(g2));	/* u dw/dx phi_n */
	d3 += *(gu + 1) * (*(ggu++)) * (*(g2));	/* v dw/dy phi_n */
	d3 += *(gu + 2) * (*(ggu++)) * (*(g2));	/* w dw/dz phi_n */

	dd1 += d1 * (*w);
	dd2 += d2 * (*w);
	dd3 += d3 * (*w);

	gu += 3;
	g2++;
	w++;
	lambda += Dim + 1;
    }

    vol = phgGeomGetVolume(u->g, e);
    values[0] = dd1 * vol;
    values[1] = dd2 * vol;
    values[2] = dd3 * vol;
    return;
}
/* **********************
 * \grad phi_j \times \phi_i
 * return 3 values
 *  */
void
phgQuadGradBasDotBas3(ELEMENT *e, DOF *s, int m, DOF *v, int n,
		      FLOAT *values, int order)
{
    int i, j, nvalues = DofTypeDim(v);
    const FLOAT *g1, *g2, *w, *lambda;
    FLOAT d1, d2, d3, dd1, dd2, dd3;
    QUAD *quad;

    assert(!SpecialDofType(s->type) && !SpecialDofType(v->type));

    if (nvalues != DofTypeDim(s))
	phgError(1, "%s:%d, dimensions mismatch: grad(%s) <==> (%s)\n",
		 __FILE__, __LINE__, s->name, v->name);

    if (order < 0)
	order = BasisOrder(s, e, m) - 1 + BasisOrder(v, e, n);
    if (order < 0)
	order = 0;
    quad = phgQuadGetQuad3D(order);

    g1 = phgQuadGetBasisGradient(e, s, m, quad);
    g2 = phgQuadGetBasisValues(e, v, n, quad);
    dd1 = 0.;
    dd2 = 0.;
    dd3 = 0.;
    lambda = quad->points;
    w = quad->weights;
    for (i = 0; i < quad->npoints; i++) {
	d1 = d2 = d3 = 0.;
	for (j = 0; j < nvalues; j++) {
	    d1 += *(g1++) * *g2;
	    d2 += *(g1++) * *g2;
	    d3 += *(g1++) * *g2;
	    g2++;
	}
	dd1 += d1 * *w;
	dd2 += d2 * *w;
	dd3 += d3 * *w;
	w++;
	lambda += Dim + 1;
    }
    values[0] = dd1 * phgGeomGetVolume(s->g, e);
    values[1] = dd2 * phgGeomGetVolume(s->g, e);
    values[2] = dd3 * phgGeomGetVolume(s->g, e);
    return;
}
Пример #5
0
/****************************************************************
 * Build RHS which is the residual of the nonlinear system.
 ***************************************************************/
void 
phgNSBuildSolverURHS(NSSolver *ns)
{
    GRID *g = ns->g;
    SIMPLEX *e;
    SOLVER *solver_u = ns->solver_u;
    int i, k, l, q, s;
    FLOAT *dt = ns->dt;
    BOOLEAN tstep_minus = (ns->u[-1] != NULL);
    VEC *vec_rhs = phgMapCreateVec(solver_u->rhs->map, 1);
    FLOAT Theta = _nsp->Theta, nu = _nsp->nu, Thet1;
    int viscosity_type = ns->viscosity_type;

    SURF_BAS *surf_bas = ns->surf_bas;
    DOF *surf_dof = surf_bas->dof;
    BOOLEAN *rotated = surf_bas->rotated;
    const FLOAT *Trans = DofData(surf_dof);


#if STEADY_STATE
    assert(fabs(Theta - 1) < 1e-12);
    Thet1 = 0; Unused(Thet1);
    Unused(dt);
#else
    Thet1 = 1 - Theta;
    Unused(dt);
#endif /* STEADY_STATE */

    phgPrintf("   DB_mask: [");
    for (k = 0; k < Dim; k++)
	phgPrintf("%d ", ns->u[1]->DB_masks[k]);
    phgPrintf("]   ");

    nu_max = -1e10;
    nu_min = +1e10;

    phgVecDisassemble(vec_rhs);
    ForAllElements(g, e) {
	int M = ns->u[1]->type->nbas;	/* num of bases of Velocity */
	int N = ns->p[1]->type->nbas;	/* num of bases of Pressure */
	int order = DofTypeOrder(ns->u[1], e) * 3 - 1; /* Note:
							*   quad order is really high here,
							*   highest order term (u \nabla u, phi)  */
	FLOAT bufu[M], bufp[N], rhsu[M][Dim], rhsp[N];
	INT Iu[M][Dim], Ip[N];
	QUAD *quad;
	FLOAT vol, area, det;
	const FLOAT *w, *p, *normal,
	    **vu, *vu_queue[3],
	    *vf[2], *gu[2], *vp[2], *vw, *vT;
	FLOAT *vf_cache[2];

	vu = vu_queue + 1;

	quad = phgQuadGetQuad3D(order);
	vu[0] = phgQuadGetDofValues(e, ns->u[0], quad);	           /* u^{n} */
	vp[0] = phgQuadGetDofValues(e, ns->p[0], quad);	           /* p^{n} */
	gu[0] = phgQuadGetDofValues(e, ns->gradu[0], quad);        /* grad u^{n} */
	vw = phgQuadGetDofValues(e, ns->wind, quad);               /* wind */
	vT = phgQuadGetDofValues(e, ns->T[1], quad);	           /* T^{n} */


	if (tstep_minus) { 
	    vu[-1] = phgQuadGetDofValues(e, ns->u[-1], quad);      /* u^{n-1} */
	} else {
	    vu[-1] = vu[0];
	}

#if STEADY_STATE || TIME_DEP_NON
	vu[1] = phgQuadGetDofValues(e, ns->u[1], quad);          /* u^{n+1} */
	gu[1] = phgQuadGetDofValues(e, ns->gradu[1], quad);      /* grad u^{n} */
	vp[1] = phgQuadGetDofValues(e, ns->p[1], quad);          /* p^{n+1} */
#else
	TIME_DEP_LINEAR_ENTRY; /* Unavailable */
#endif /* STEADY_STATE || TIME_DEP_NON */

	Unused(l);
	Unused(vf); Unused(vf_cache); 

	if (!_nsp->no_extern_source) {
	    /* cache f values */
	    for (l = 0; l < 2; l++) {
		const FLOAT *cache;
		size_t cache_size;
		setFuncTime(ns->time[l]); /* set static time in ins-test.c */

		/* cache f */
		cache_size = Dim * quad->npoints * sizeof(FLOAT);
		cache = phgQuadGetFuncValues(g, e, Dim, func_f, quad);
		vf[l] = vf_cache[l] = phgAlloc(cache_size);
		memcpy(vf_cache[l], cache, cache_size);

		phgQuadGetFuncValues(NULL, NULL, 0, NULL, NULL); /* clear cache */
	    }
	}

	/* Global Matrix */
	Bzero(rhsu); Bzero(rhsp);
    
	p = quad->points;
	w = quad->weights;
	for (q = 0; q < quad->npoints; q++) {
	    phgGeomGetCurvedJacobianAtLambda(g, e, p, &det);
	    vol = fabs(det / 6.);

	    /* rhs u */
	    for (i = 0; i < M; i++) {
		/* interior node or Neumann */
		const FLOAT *gi_u = phgQuadGetBasisValues(e, ns->u[1], i, quad) + q;       /* phi_i */
		const FLOAT *ggi_u = phgQuadGetBasisCurvedGradient(e, ns->u[1], i, quad, q);    /* grad phi_i */

		for (k = 0; k < Dim; k++) {
#if ICE_BENCH_TEST
		    nu = get_effective_viscosity(gu[1], 0, 0, viscosity_type);
		    FLOAT eu[DDim];

		    MAT3_SYM(gu[1], eu);
		    rhsu[i][k] += vol*(*w) * EQU_SCALING * (- nu * INNER_PRODUCT(eu+k*Dim, ggi_u) 
							    + (*vp[1]) * *(ggi_u+k) * LEN_SCALING * PRES_SCALING
							    );     /* left */

		    if (k == Z_DIR) { 
			const FLOAT rho = RHO_ICE;
			const FLOAT grav = GRAVITY;
			const FLOAT a = SEC_PER_YEAR;
			const FLOAT f = rho*grav * EQU_SCALING * LEN_SCALING2; 

			Unused(a);
			rhsu[i][k] += vol*(*w) * (-f * (*gi_u) 
						  ); /* right */
		    }


#elif ESIMINT_TEST ||				\
    HEINO_TEST ||				\
    TEST_CASE == ICE_GREEN_LAND
		    nu = get_effective_viscosity(gu[1], *vT, 0, viscosity_type);
		    FLOAT eu[DDim];

		    MAT3_SYM(gu[1], eu);
		    rhsu[i][k] += vol*(*w) * EQU_SCALING * (- nu * INNER_PRODUCT(eu+k*Dim, ggi_u) 
							    + (*vp[1]) * *(ggi_u+k) * LEN_SCALING * PRES_SCALING
							    );     /* left */

		    if (k == Z_DIR) { 
			const FLOAT rho = RHO_ICE;
			const FLOAT grav = GRAVITY;
			const FLOAT a = SEC_PER_YEAR;
			const FLOAT f = rho*grav * EQU_SCALING * LEN_SCALING2; 

			Unused(a);
			rhsu[i][k] += vol*(*w) * (-f * (*gi_u) 
						  ); /* right */
		    }

#elif STEADY_STATE
		    rhsu[i][k] += vol*(*w) * (- nu * INNER_PRODUCT(gu[1]+k*Dim, ggi_u)
					      + (*vp[1]) * *(ggi_u+k)
					      );     /* left */
		    if (!_nsp->no_extern_source)
			rhsu[i][k] += vol*(*w) * (*(vf[1]+k) * (*gi_u)
						  ); /* right */
#elif TIME_DEP_NON
		    rhsu[i][k] -= vol*(*w) * ((vu[1][k] - vu[0][k]) * (*gi_u) / dt[0]
					      + Theta * (nu * INNER_PRODUCT(gu[1]+k*Dim, ggi_u)
							 )
					      - (*vp[1]) * *(ggi_u+k)
					      + Thet1 * (nu * INNER_PRODUCT(gu[0]+k*Dim, ggi_u)
							 )
					      );     /* left */
		    if (!_nsp->no_extern_source)
			rhsu[i][k] += vol*(*w) * (Theta * *(vf[1]+k) * (*gi_u)
						  + Thet1 * *(vf[0]+k) * (*gi_u)
						  ); /* right */
#else
		    TIME_DEP_LINEAR_ENTRY; /* Unavailable */
#endif /* STEADY_STATE */
		}
	    }

	    /* rhs p */
	    for (i = 0; i < N; i++) {
		const FLOAT *gi_p = phgQuadGetBasisValues(e, ns->p[1], i, quad) + q;       /* psi_i */
		FLOAT divu1 = gu[1][0] + gu[1][4] + gu[1][8];
		//FLOAT divu0 = gu[0][0] + gu[0][4] + gu[0][8];
		rhsp[i] += vol*(*w) * (divu1 * (*gi_p)
				       );
	    }
	    
	    if (tstep_minus) 
		vu[-1] += Dim;

#if STEADY_STATE || TIME_DEP_NON
	    vu[1] += Dim;
	    gu[1] += Dim*Dim;
	    vp[1]++;
#else
	    TIME_DEP_LINEAR; /* Unavailable */
#endif /* STEADY_STATE || TIME_DEP_NON */
	    vu[0] += Dim;
	    gu[0] += Dim * Dim;
	    vp[0]++; 
	    vw += Dim;
	    if (!_nsp->no_extern_source) {
		vf[0] += Dim; vf[1] += Dim;
	    }
	    vT++;
	    w++; p += Dim + 1;
	}

	if (!_nsp->no_extern_source) {
	    phgFree(vf_cache[0]);
	    phgFree(vf_cache[1]);
	}

	normal = NULL; Unused(normal);
	area = 0; Unused(area);

	if (!_nsp->enclosed_flow) {
	    /* slip boundary */
	    for (s = 0; s < NFace; s++) {
		if (e->bound_type[s] & INFLOW) {
		    int v0, v1, v2;
		    int nbas_face = NbasFace(ns->u[1]);
		    SHORT bases[nbas_face];
		    FLOAT lambda[Dim + 1], x,y,z, beta;
		    order = DofTypeOrder(ns->u[1], e) * 3 - 1;

		    phgDofGetBasesOnFace(ns->u[1], e, s, bases);
		    v0 = GetFaceVertex(s, 0);
		    v1 = GetFaceVertex(s, 1);
		    v2 = GetFaceVertex(s, 2);
		    lambda[s] = 0.;
		    
		    area = phgGeomGetFaceArea(g, e, s);
		    normal = phgGeomGetFaceOutNormal(g, e, s);
		    quad = phgQuadGetQuad2D(order);
		
		    p = quad->points;
		    w = quad->weights;
		    for (q = 0; q < quad->npoints; q++) {
			FLOAT vu[Dim];
			lambda[v0] = *(p++);
			lambda[v1] = *(p++);
			lambda[v2] = *(p++);
			
			phgGeomLambda2XYZ(g, e, lambda, &x, &y, &z);
			func_beta(x, y, z, &beta);
			phgDofEval(ns->u[1], e, lambda, vu);

			for (i = 0; i < nbas_face; i++) {
			    int ii = bases[i];
			    FLOAT gi_u = 
				*ns->u[1]->type->BasFuncs(ns->u[1], e, ii, ii + 1, lambda);

			    for (k = 0; k < Dim; k++) {
#if STEADY_STATE
				rhus[ii][k] += 0.;
#elif TIME_DEP_NON
#  if USE_SLIDING_BC
				abort();
				rhsu[ii][k] += SIGN_FRICTION * area*(*w) * beta * vu[k] * (gi_u)
				    * EQU_SCALING * LEN_SCALING;
#  else
				Unused(gi_u);
#  endif
#else
				TIME_DEP_LINEAR_ENTRY; /* Unavailable */
#endif /* STEADY_STATE */
			    }
			}     /* end of bas_i */
			w++;
		    }		/* end of quad point */
		}		/* end of face outflow */
	    }			/* end of all outflow face in element */
	}                       /* end out flow boundary */


#if USE_SLIDING_BC
	/* Rotate bases */
	for (i = 0; i < M; i++) {
	    INT id = phgDofMapE2D(surf_dof, e, i * (Dim*Dim)) / (Dim*Dim);
	    if (!rotated[id])
		continue;	
	    const FLOAT *trans = Trans + id*(Dim*Dim);

	    trans_left(&rhsu[i][0], 1, 1, trans);
	}
#else
	Unused(Trans);
	Unused(rotated);
#endif


	/* Map: Element -> system */
	for (i = 0; i < M; i++)
	    for (k = 0; k < Dim; k++)
		Iu[i][k] = phgMapE2L(solver_u->rhs->map, 0, e, i * Dim + k);
	for (i = 0; i < N; i++)
	    Ip[i] = phgMapE2L(solver_u->rhs->map, 1, e, i);

	/* set velocity dirichlet bdry */
	FLOAT tmp[Dim];
	for (i = 0; i < M; i++)
	    for (k = 0; k < Dim; k++)
		if (phgDofDirichletBC_(ns->u[1], e, i*Dim+k, NULL, bufu, tmp,
				       DOF_PROJ_NONE)) {
		    rhsu[i][k] = 0.;
		}

#if STEADY_STATE || TIME_DEP_NON
	/* set pressure dirichlet bdry for pinned point */
	for (i = 0; i < N; i++)
	    if (phgDofDirichletBC(ns->p[1], e, i, NULL, bufp, &rhsp[i],
				  DOF_PROJ_NONE)) {
		if (!_nsp->enclosed_flow)
		    phgError(1, "No dirichlet bc for Unenclosed flow!\n");
		if (_nsp->pin_node) {
# if PIN_AT_ROOT 
		    if (g->rank != 0)
		    	phgError(1, "Pinned node only on rank 0!\n");
		    if (g, e->verts[i] != ns->pinned_node_id)
			phgError(1, "Build rhs: pinned node e:%d, bas:%d, [%d] and [%d] "
				 "doesn't coincide when build RHS!\n", 
				 e->index, i, e->verts[i], ns->pinned_node_id);
# else
		    if (GlobalVertex(g, e->verts[i]) != ns->pinned_node_id)
			phgError(1, "Build rhs: pinned node e:%d, bas:%d, [%d] and [%d] "
				 "doesn't coincide when build RHS!\n", 
				 e->index, i, e->verts[i], ns->pinned_node_id);
# endif /* PIN_AT_ROOT */
		}
	    } 
#else
	TIME_DEP_LINEAR; /* Unavailable */
#endif /* STEADY_STATE || TIME_DEP_NON */


	/* Global res */
	phgVecAddEntries(vec_rhs, 0, M * Dim, Iu[0], &rhsu[0][0]);
	phgVecAddEntries(vec_rhs, 0, N, Ip, rhsp);
    }				/* end element */
Пример #6
0
void
phgNSBuildPc(NSSolver *ns)
{
    GRID *g = ns->g;
    SIMPLEX *e;
    FLOAT *dt = ns->dt;
    int i, j, q, s, k, l;
    FLOAT Theta = _nsp->Theta, nu = _nsp->nu, Thet1, nu0 = 0;
    DOF *tmp_u1 = phgDofNew(g, _nsp->utype, Dim, "tmp u1", func_u);
    int viscosity_type = ns->viscosity_type;
    LTYPE ltype = ns->ltype;


#if STEADY_STATE
    assert(fabs(Theta - 1) < 1e-12);
    Thet1 = 0; Unused(Thet1);
    Unused(dt);
#else
    Thet1 = 1 - Theta;
#endif /* STEADY_STATE */


    ForAllElements(g, e) {
	int M = ns->u[1]->type->nbas;	/* num of bases of Velocity */
	int N = ns->p[1]->type->nbas;	/* num of bases of Pressure */
	int order = 2 * DofTypeOrder(ns->p[1], e) + 
	    DofTypeOrder(ns->u[1], e) - 1; 	/* highest order term (u \nabla p, psi)  */
	FLOAT Ap[N][N], Fp[N][N], Qp[N][N], bufp[N], rhs1 = 1;
	FLOAT F[M*Dim][M*Dim], B[N][M*Dim], Bt[M*Dim][N];
	INT Ip[N];
	QUAD *quad;
	FLOAT vol, det;
	const FLOAT *w, *p, *vw, *gu, *vTe;

	quad = phgQuadGetQuad3D(order);
	vw = phgQuadGetDofValues(e, ns->wind, quad);  /* value wind */
	gu = phgQuadGetDofValues(e, ns->gradu[1], quad);        /* grad u^{n+1} */
	if (ns_params->noniter_temp)
	    vTe = phgQuadGetDofValues(e, ns->T[1], quad);  /* value temp */
	else
	    vTe = phgQuadGetDofValues(e, ns->T[0], quad);  /* value temp */
	
	vol = 0;
	Bzero(Ap); Bzero(Fp); Bzero(Qp); 
	Bzero(F); Bzero(Bt); Bzero(B);
	Bzero(bufp); 

	p = quad->points;
	w = quad->weights;
	for (q = 0; q < quad->npoints; q++) {
	    phgGeomGetCurvedJacobianAtLambda(g, e, p, &det);
	    vol = fabs(det / 6.);

	    for (i = 0; i < N; i++) {
		const FLOAT *gi = phgQuadGetBasisValues(e, ns->p[1], i, quad) + q;       /* phi_i */
		const FLOAT *ggi = phgQuadGetBasisCurvedGradient(e, ns->p[1], i, quad, q);    /* grad phi_i */
		for (j = 0; j < N; j++) {
		    const FLOAT *gj = phgQuadGetBasisValues(e, ns->p[1], j, quad) + q;       /* phi_j */
		    const FLOAT *ggj = phgQuadGetBasisCurvedGradient(e, ns->p[1], j, quad, q);    /* grad phi_i */
		    
		    nu = get_effective_viscosity(gu, *vTe, 0, viscosity_type);
		    if (i == 0 && j == 0)
			nu0 += nu;
#if ICE_BENCH_TEST ||				\
    ESIMINT_TEST ||				\
    HEINO_TEST ||				\
    TEST_CASE == ICE_EXACT	||		\
    TEST_CASE == ICE_GREEN_LAND
		    Unused(dt);
		    /* Note: B Q^-1 Bt ~ Ap(nu=1),
		     *       Fp(nu varies) is very different to Ap */
		    Ap[i][j] += vol*(*w) * INNER_PRODUCT(ggj, ggi);
#  if USE_QP_ONLY

		    //Qp[i][j] += vol*(*w) * LEN_SCALING * PRES_SCALING /(nu) * (*gj) * (*gi);
		    Qp[i][j] += vol*(*w) * 1. /(EQU_SCALING * nu) * (*gj) * (*gi);
		    /* if (i < NVert && j < NVert) { */
		    /* 	Qp[i][j] += vol*(*w) * LEN_SCALING * PRES_SCALING / (nu) * (*gj) * (*gi); */
		    /* } else if (i == NVert && j == NVert) { */
		    /* 	Qp[i][j] += vol*(*w) * LEN_SCALING * PRES_SCALING / (nu) * (*gj) * (*gi); */
		    /* } */

#  else
		    Qp[i][j] += vol*(*w) * (*gj) * (*gi);
#  endif
		    Fp[i][j] += vol*(*w) * (EQU_SCALING * nu * INNER_PRODUCT(ggj, ggi)
					    );
#elif STEADY_STATE 
		    Ap[i][j] += vol*(*w) * INNER_PRODUCT(ggj, ggi);
		    Qp[i][j] += vol*(*w) * (*gj) * (*gi);
		    Fp[i][j] += vol*(*w) * (nu * INNER_PRODUCT(ggj, ggi) * EQU_SCALING
					    );
#elif TIME_DEP_NON
		    Ap[i][j] += vol*(*w) * INNER_PRODUCT(ggj, ggi);
		    Qp[i][j] += vol*(*w) * (*gj) * (*gi);
		    Fp[i][j] += vol*(*w) * ((*gj) * (*gi) / dt[0]
					    + Theta * (nu * INNER_PRODUCT(ggj, ggi)
						       )
					    );
#else
		    TIME_DEP_LINEAR_ENTRY; /* Unavailable */
#endif /* STEADY_STATE */
		}
	    }

	    vw += Dim; 
	    gu += DDim;
	    vTe++;
	    w++; p += Dim+1;
	}


	/* Map: Element -> system */
	for (i = 0; i < N; i++) 
	    Ip[i] = phgMapE2L(_pcd->matFp->cmap, 0, e, i);

	/*
	 * PCD boundary setup I:
	 * Automaticly decide inflow boundary condition using wind direction.
	 *
	 * NOT active.
	 * */
	if (FALSE && !_nsp->pin_node) {
	    for (i = 0; i < N; i++) {
		BOOLEAN flag_inflow = FALSE;
		for (s = 0; s < NFace; s++) {
		    SHORT bases[NbasFace(ns->p[1])];
		    FLOAT *coord, vw_[3]; 
		    const FLOAT *lam, *normal;

		    if (!(e->bound_type[s] & BDRY_MASK))
			//if (!(e->bound_type[s] & INFLOW))
			continue;	/* boundary face */

		    phgDofGetBasesOnFace(ns->p[1], e, s, bases);
		    for (j = 0; j < NbasFace(ns->p[1]); j++) 
			if (i == bases[j]) {
			    normal = phgGeomGetFaceOutNormal(g, e, s);
			    coord = phgDofGetElementCoordinates(ns->p[1], e, i);
			    lam = phgGeomXYZ2Lambda(g, e, coord[0], coord[1], coord[2]);
			    phgDofEval(tmp_u1, e, lam, vw_);
			    if (INNER_PRODUCT(vw_, normal) > 1e-8) 
				flag_inflow = TRUE;
			}
		}
		
		if (flag_inflow) {
		    Bzero(bufp); bufp[i] = 1.0;
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp);
		    //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp);
		    phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1);
		}
		else {
		    /* interior node Or Neumann */
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, Ap[i]);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]);
		    //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]);
		}
		phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]);
	    }
	} 
	/*
	 * PCD boundary setup II:
	 * Enclose flow: use pinnode boundary.
	 *
	 * Qp is pinned, this is different to open flow.
	 * 
	 * */
	else if (_nsp->pin_node) {

	    for (i = 0; i < N; i++) {
		if (phgDofDirichletBC(_pcd->pbc, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)) {
#if PIN_AT_ROOT 
		    if (g->rank != 0)
		    	phgError(1, "Pinned node only on rank 0!\n");
		    if (e->verts[i] != ns->pinned_node_id)
			phgError(1, "pinned node [%d] & [%d] doesn't coincide when build pc!\n",
				 e->verts[i], ns->pinned_node_id);
#else
		    if (GlobalVertex(g, e->verts[i]) != ns->pinned_node_id)
			phgError(1, "pinned node [%d] & [%d] doesn't coincide when build pc!\n",
				 e->verts[i], ns->pinned_node_id);
#endif /* PIN_AT_ROOT */

		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp);
		    phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1);
		} else {
		    /* interior node Or Neumann */
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, Ap[i]);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]);
		    phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]);
		}
	    }
	}
	/*
	 * PCD boundary setup III:
	 * Open flow: there could be varies kinds of combination on seting up
	 *   boundary conditon, but Inflow:Robin & Outflow:scaled Dirich is
	 *   prefered. See Ref[2].
	 * 
	 * */
	else {
	    for (i = 0; i < N; i++) {

		/*****************/
                /* Inflow	 */
                /*****************/
#warning PCD B.C.: Step 2.1. build mat, all neumann, add dirich entries
		if (FALSE && phgDofDirichletBC(_pcd->dof_inflow, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)) {
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp);
		    phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1);
		} else if (FALSE && phgDofDirichletBC(_pcd->dof_outflow, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)
			   && !(phgDofGetElementBoundaryType(ns->p[1], e, i) & INFLOW) ) {

		    ERROR_MSG("Fp, Qp");
		    nu = get_effective_viscosity(NULL, 0, 0, viscosity_type);
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp);
		    bufp[i] *= EQU_SCALING * nu;
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp);
		    phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1);

		    //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp);
		} else if (FALSE && phgDofDirichletBC(_pcd->pbc, e, i, NULL, bufp, NULL, DOF_PROJ_NONE)) {
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, bufp);
		    phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, bufp);
		    phgVecAddEntries(_pcd->rhsScale, 0, 1, Ip + i, &rhs1);
		}
		else if (FALSE) {
		    /* interior node Or Neumann */

		    ERROR_MSG("Fp, Qp");
		    phgMatAddEntries(_pcd->matAp, 1, Ip + i, N, Ip, Ap[i]);
		    phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]);
		    //phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]);
		}

		/******************/
                /* No bdry	  */
                /******************/
		//phgMatAddEntries(_pcd->matFp, 1, Ip + i, N, Ip, Fp[i]);
		phgMatAddEntries(_pcd->matQp, 1, Ip + i, N, Ip, Qp[i]);
	    }
	}


	if (0) {
	    /* Special term <[[p_i]], [[p_j]]> */
	    int face;
	    nu0 /= quad->npoints;
	    for (face = 0; face < NFace; face++) {
		FLOAT area =  phgGeomGetFaceArea(g, e, face);
		//FLOAT value = {area, -area};
		FLOAT values[2] = {vol * 1. /(EQU_SCALING * nu0),
				   -vol * 1. /(EQU_SCALING * nu0)};
		SIMPLEX *e_neigh;

		phgMatAddEntries(_pcd->matQp, 1, Ip+NVert, 1, Ip+NVert, values);
		if ((e_neigh = GetNeighbour(e, face)) != NULL) {
		    INT Ip_neigh = phgMapE2L(_pcd->matFp->cmap, 0, e_neigh, NVert);
		    phgMatAddEntries(_pcd->matQp, 1, Ip+NVert, 1, &Ip_neigh, values + 1);
		}
	    }
	}

    }	/* end element */