/* Scatter the points from u to U using the map in B->bmap */ static void ScatrBndry_Stokes(double *u, Element_List **V, Bsystem **B) { register int j,k; int eDIM = V[0]->flist[0]->dim(); int nel = B[0]->nel, Nbmodes, *bmap; double *hj; double *sc = B[0]->signchange; for(k = 0; k < nel; ++k) { bmap = B[eDIM]->bmap[k]; Nbmodes = V[0]->flist[k]->Nbmodes; for(j = 0; j < eDIM; ++j) { hj = V[j]->flist[k]->vert->hj; dgathr(Nbmodes, u, bmap + j*Nbmodes, hj); dvmul (Nbmodes, sc, 1, hj, 1, hj,1); } V[eDIM]->flist[k]->vert[0].hj[0] = u[bmap[eDIM*Nbmodes]]; sc += Nbmodes; } }
void Quad::fillvec(Mode *v, double *f){ register int i; for(i = 0; i < qb; ++i) dcopy(qa,v->a,1,f+i*qa,1); for(i = 0; i < qa; ++i) dvmul(qb,v->b,1,f+i,qa,f+i,qa); }
static void solve_pressure(Element_List **V, double *rhs, Bsystem **Vbsys) { register int i; int eDIM = V[0]->flist[0]->dim(); int info,N,nbl,nblv,id,*bmap; Bsystem *B = Vbsys[0], *PB = Vbsys[eDIM]; Element *E; double *hj,*tmp; double *sc = B->signchange; if(eDIM == 2) tmp = dvector(0,8*LGmax); else tmp = dvector(0,18*(LGmax-1)*(LGmax-1)); /* back solve for pressure */ for(E=V[eDIM]->fhead;E;E=E->next) { N = E->Nmodes - 1; id = E->geom->id; hj = E->vert->hj + 1; E->state = 't'; /* solve interior and negate */ if(PB->lambda->wave) { dgetrs('N', N, 1, PB->Gmat->invc[id], N, PB->Gmat->cipiv[id],hj, N, info); } else { dpptrs('L', N, 1, PB->Gmat->invc[id], hj, N, info); dneg(N,hj,1); } bmap = PB->bmap[E->id]; nblv = V[0]->flist[E->id]->Nbmodes; nbl = eDIM*nblv+1; for(i = 0; i < nbl; ++i) tmp[i] = rhs[bmap[i]]; for(i = 0; i < eDIM; ++i) dvmul(nblv,sc,1,tmp+nblv*i,1,tmp+nblv*i,1); if(PB->lambda->wave) dgemv('N', N, nbl,-1.0, PB->Gmat->invcd[id], N, tmp, 1, 1.,hj,1); else dgemv('N', N, nbl,-1.0, PB->Gmat->binvc[id], N, tmp, 1, 1.,hj,1); sc += nblv; } free(tmp); }
void A_Stokes(Element_List **V, Bsystem **B, double *p, double *w, double **wk) { int eDIM = V[0]->flist[0]->dim(); int nel = B[0]->nel, Nbmodes; double **a = B[eDIM]->Gmat->a; int **bmap = B[eDIM]->bmap; double *sc = B[0]->signchange; register int i,j,k; dzero(B[eDIM]->nsolve,w,1); for(k = 0; k < nel; ++k) { Nbmodes = V[0]->flist[k]->Nbmodes; /* Put p boundary modes into wk[0] and impose continuity */ dgathr(eDIM*Nbmodes+1, p, bmap[k], wk[0]); for(j = 0; j < eDIM; ++j) dvmul (Nbmodes, sc, 1, wk[0] + j*Nbmodes, 1, wk[0] + j*Nbmodes,1); dspmv('U',eDIM*Nbmodes+1,1.0,a[V[0]->flist[k]->geom->id], wk[0],1,0.0,wk[1],1); /* do sign change back and put into w */ for(j = 0; j < eDIM; ++j) for(i = 0; i < Nbmodes; ++i) w[bmap[k][i+j*Nbmodes]] += sc[i]*wk[1][i+j*Nbmodes]; w[bmap[k][eDIM*Nbmodes]] += wk[1][eDIM*Nbmodes]; sc += Nbmodes; } /* set first pressure constant to zero */ if(B[eDIM]->singular) w[eDIM*B[0]->nsolve] = 0.0; }
double Quad::get_1diag_massmat(int ID){ double *wa, *wb; double *wvec = dvector(0, qtot-1), vmmat; Mode mw,*m; #ifndef PCONTBASE double **ba, **bb; Mode m1; get_moda_GL (qa, &ba); get_moda_GL (qb, &bb); m1.a = ba[ID]; m1.b = bb[ID]; m = &m1; #else Basis *b = getbasis(); m = b->vert+ID; #endif getzw(qa,&wa,&wa,'a'); getzw(qb,&wb,&wb,'a'); mw.a = wa; mw.b = wb; fillvec(&mw, wvec); if(curvX) dvmul(qtot, wvec, 1, geom->jac.p, 1, wvec, 1); else dscal(qtot, geom->jac.d, wvec, 1); vmmat = Quad_mass_mprod(this, m, wvec); free(wvec); return vmmat; }
// ============================================================================ void Tri::BET_Mat(Element *P, LocMatDiv *bet, double *beta, double *sigma) { register int i,j,n; const int nbl = Nbmodes, N = Nmodes - Nbmodes; int L; Basis *b = getbasis(); double **dxb = bet->Dxb, // MSB: dx corresponds to bar(beta) **dxi = bet->Dxi, **dyb = bet->Dyb, // MSB: dy corresponds to sigma **dyi = bet->Dyi; char orig_state = state; /* fill boundary systems */ for(i = 0,n=0; i < Nverts; ++i,++n) { fillElmt(b->vert+i); dvmul(qtot,beta,1,*h,1,*P->h,1); #ifndef PCONTBASE P->Ofwd(*P->h,P->vert->hj,P->dgL); #else P->Iprod(P); #endif dcopy(P->Nmodes,P->vert->hj,1,*dxb + n,nbl); dvmul(qtot,sigma,1,*h,1,*P->h,1); #ifndef PCONTBASE P->Ofwd(*P->h,P->vert->hj,P->dgL); #else P->Iprod(P); #endif dcopy(P->Nmodes,P->vert->hj,1,*dyb + n,nbl); } for(i = 0; i < Nedges; ++i) for(j = 0; j < edge[i].l; ++j, ++n) { fillElmt(b->edge[i]+j); dvmul(qtot,beta,1,*h,1,*P->h,1); #ifndef PCONTBASE P->Ofwd(*P->h,P->vert->hj,P->dgL); #else P->Iprod(P); #endif dcopy(P->Nmodes,P->vert->hj,1,*dxb + n,nbl); dvmul(qtot,sigma,1,*h,1,*P->h,1); #ifndef PCONTBASE P->Ofwd(*P->h,P->vert->hj,P->dgL); #else P->Iprod(P); #endif dcopy(P->Nmodes,P->vert->hj,1,*dyb + n,nbl); } L = face->l; for(i = 0,n=0; i < L;++i) for(j = 0; j < L-i; ++j,++n) { fillElmt(b->face[0][i]+j); dvmul(qtot,beta,1,*h,1,*P->h,1); #ifndef PCONTBASE P->Ofwd(*P->h,P->vert->hj,P->dgL); #else P->Iprod(P); #endif dcopy(P->Nmodes,P->vert->hj,1,*dxi + n,N); dvmul(qtot,sigma,1,*h,1,*P->h,1); #ifndef PCONTBASE P->Ofwd(*P->h,P->vert->hj,P->dgL); #else P->Iprod(P); #endif dcopy(P->Nmodes,P->vert->hj,1,*dyi + n,N); } state = orig_state; /* negate all systems to that the whole operator can be treated as positive when condensing */ /* dneg(nbl*P->Nmodes,*dxb,1); dneg(nbl*P->Nmodes,*dyb,1); dneg(N *P->Nmodes,*dxi,1); dneg(N *P->Nmodes,*dyi,1); */ }
// ============================================================================ void Tri::PSE_Mat(Element *E, LocMat *pse, double *DU) { register int i,j,n; int nbl, L, N, Nm, qt, asize = pse->asize, csize = pse->csize; Basis *B = E->getbasis(); double *Eh; double **a = pse->a; double **b = pse->b; double **c = pse->c; double **d = pse->d; // This routine is placed within the for loop for each element // at around line 30 in the code. Therefore, this is within an element nbl = E->Nbmodes; N = E->Nmodes - nbl; Nm = E->Nmodes; qt = E->qa*E->qb; Eh = dvector(0, qt - 1); // Temporary storage for E->h[0] ------ dcopy(qt, E->h[0], 1, Eh, 1); // Fill A and D ---------------------------------------------------- for(i = 0, n = 0; i < E->Nverts; ++i, ++n) { E->fillElmt(B->vert + i); // ROUTINE THAT DOES INNERPRODUCT AND EVALUATES TERM-------------- dvmul(qt, DU, 1, E->h[0], 1, E->h[0], 1); E->Iprod(E); // --------------------------------------------------------------- dcopy(nbl, E->vert->hj, 1, *a + n, asize); dcopy(N, E->vert->hj + nbl, 1, d[n], 1); } dcopy(qt, Eh, 1, E->h[0], 1); for(i = 0; i < E->Nedges; ++i) { for(j = 0; j < (L = E->edge[i].l); ++j, ++n) { E->fillElmt(B->edge[i] + j); // ROUTINE THAT DOES INNERPRODUCT AND EVALUATES TERM------------ dvmul(qt, DU, 1, E->h[0], 1, E->h[0], 1); E->Iprod(E); // ------------------------------------------------------------- dcopy(nbl, E->vert->hj, 1, *a + n, asize); dcopy(N, E->vert->hj + nbl, 1, d[n], 1); } } dcopy(qt, Eh, 1, E->h[0], 1); // Fill B and C ---------------------------------------------------- L = E->face->l; for(i = 0, n = 0; i < L; ++i) { for(j = 0; j < L-i; ++j, ++n) { E->fillElmt(B->face[0][i]+j); // ROUTINE THAT DOES INNERPRODUCT AND EVALUATES TERM----------- dvmul(qt, DU, 1, E->h[0], 1, E->h[0], 1); E->Iprod(E); // ------------------------------------------------------------ dcopy(nbl, E->vert->hj , 1, *b + n, csize); dcopy(N, E->vert->hj + nbl, 1, *c + n, csize); } } // ----------------------------------------------------------------- free(Eh); }
void Tri_setup_iprodlap() { int qa = QGmax; int qb; if(LZero) // use extra points if have Legendre zeros in 'b' direction qb = QGmax; else qb = QGmax-1; int qc = 0; int L = LGmax; // Set up dummy element with maximum quadrature/edge order Coord X; X.x = dvector(0,NTri_verts-1); X.y = dvector(0,NTri_verts-1); X.x[0] = 0.0; X.x[1] = 1.0; X.x[2] = 0.0; X.y[0] = 0.0; X.y[1] = 0.0; X.y[2] = 1.0; Tri *T = (Tri*) new Tri(0,'T', L, qa, qb, qc, &X); free(X.x); free(X.y); int i,j,k,n; int facs = Tri_DIM*Tri_DIM; Tri_ipmodes = T->Nmodes; Tri_iplap = (double***) calloc(facs,sizeof(double**)); for(i = 0; i < facs; ++i) { Tri_iplap[i] = dmatrix(0, T->Nmodes-1, 0, T->Nmodes-1); dzero(T->Nmodes*T->Nmodes, Tri_iplap[i][0], 1); } // Set up gradient basis Basis *B,*DB; Mode *w,*m,*m1,*md,*md1,**gb,**gb1,*fac; double *z; B = T->getbasis(); DB = T->derbasis(); fac = B->vert; w = mvector(0,0); gb = (Mode **) malloc(T->Nmodes*sizeof(Mode *)); gb[0] = mvecset(0,Tri_DIM*T->Nmodes,qa, qb, qc); gb1 = (Mode **) malloc(T->Nmodes*sizeof(Mode *)); gb1[0] = mvecset(0,Tri_DIM*T->Nmodes,qa, qb, qc); for(i = 1; i < T->Nmodes; ++i) gb[i] = gb[i-1]+Tri_DIM; for(i = 1; i < T->Nmodes; ++i) gb1[i] = gb1[i-1]+Tri_DIM; getzw(qa,&z,&w[0].a,'a'); getzw(qb,&z,&w[0].b,'b'); /* fill gb with basis info for laplacian calculation */ // vertex modes m = B->vert; md = DB->vert; for(i = 0,n=0; i < T->Nverts; ++i,++n) T->fill_gradbase(gb[n],m+i,md+i,fac); // edge modes for(i = 0; i < T->Nedges; ++i) { m1 = B ->edge[i]; md1 = DB->edge[i]; for(j = 0; j < T->edge[i].l; ++j,++n) T->fill_gradbase(gb[n],m1+j,md1+j,fac); } // face modes for(i = 0; i < T->Nfaces; ++i) for(j = 0; j < T->face[0].l; ++j) { m1 = B ->face[i][j]; md1 = DB->face[i][j]; for(k = 0; k < T->face[0].l-j; ++k,++n) T->fill_gradbase(gb[n],m1+k,md1+k,fac); } /* multiply by weights */ for(i = 0; i < T->Nmodes; ++i) { Tri_mvmul2d(qa,qb,qc,gb[i] ,w,gb1[i]); Tri_mvmul2d(qa,qb,qc,gb[i]+1,w,gb1[i]+1); } // Calculate Laplacian inner products double s1, s2, s3, s4, s5, s6, s7, s8; double *tmp = dvector(0, QGmax-1); fac = B->vert+1; for(i = 0; i < T->Nmodes; ++i) for(j = 0; j < T->Nmodes; ++j) { s1 = ddot(qa,gb[i][0].a,1,gb1[j][0].a,1); s1 *= ddot(qb,gb[i][0].b,1,gb1[j][0].b,1); dvmul(qa, fac->a, 1, gb[i][0].a, 1, tmp, 1); s2 = ddot(qa, tmp,1,gb1[j][0].a,1); s2 *= ddot(qb,gb[i][0].b,1,gb1[j][0].b,1); s3 = ddot(qa, tmp,1,gb1[j][1].a,1); s3 *= ddot(qb,gb[i][0].b,1,gb1[j][1].b,1); dvmul(qa, fac->a, 1, tmp, 1, tmp, 1); s4 = ddot(qa, tmp,1,gb1[j][0].a,1); s4 *= ddot(qb,gb[i][0].b,1,gb1[j][0].b,1); s5 = ddot(qa,gb[i][0].a,1,gb1[j][1].a,1); s5 *= ddot(qb,gb[i][0].b,1,gb1[j][1].b,1); s6 = ddot(qa,gb[i][1].a,1,gb1[j][0].a,1); s6 *= ddot(qb,gb[i][1].b,1,gb1[j][0].b,1); dvmul(qa, fac->a, 1, gb[i][1].a, 1, tmp, 1); s7 = ddot(qa, tmp,1,gb1[j][0].a,1); s7 *= ddot(qb,gb[i][1].b,1,gb1[j][0].b,1); s8 = ddot(qa,gb[i][1].a,1,gb1[j][1].a,1); s8 *= ddot(qb,gb[i][1].b,1,gb1[j][1].b,1); Tri_iplap[0][i][j] = s1; Tri_iplap[1][i][j] = s5 + 2*s2 + s6; Tri_iplap[2][i][j] = s3 + s4 + s7 + s8; } /* fill gb with basis info for mass matrix calculation */ // vertex modes m = B->vert; for(i = 0,n=0; i < T->Nverts; ++i,++n) { dcopy(qa, m[i].a, 1, gb[n]->a, 1); dcopy(qb, m[i].b, 1, gb[n]->b, 1); } // edge modes for(i = 0; i < T->Nedges; ++i) { m1 = B ->edge[i]; for(j = 0; j < T->edge[i].l; ++j,++n) { dcopy(qa, m1[j].a, 1, gb[n]->a, 1); dcopy(qb, m1[j].b, 1, gb[n]->b, 1); } } // face modes for(i = 0; i < T->Nfaces; ++i) for(j = 0; j < T->face[i].l; ++j) { m1 = B ->face[i][j]; for(k = 0; k < T->face[i].l-j; ++k,++n) { dcopy(qa, m1[k].a, 1, gb[n]->a, 1); dcopy(qb, m1[k].b, 1, gb[n]->b, 1); } } /* multiply by weights */ for(i = 0; i < T->Nmodes; ++i) Tri_mvmul2d(qa,qb,qc,gb[i] ,w,gb1[i]); for(i = 0; i < T->Nmodes; ++i) for(j = 0; j < T->Nmodes; ++j) { s1 = ddot(qa,gb[i][0].a,1,gb1[j][0].a,1); s1 *= ddot(qb,gb[i][0].b,1,gb1[j][0].b,1); Tri_iplap[3][i][j] = s1; } free_mvec(gb[0]) ; free((char *) gb); free_mvec(gb1[0]); free((char *) gb1); free((char*)w); free(tmp); delete (T); }
void Precon_Stokes(Element_List *U, Bsystem *B, double *r, double *z){ switch(B->Precon){ case Pre_Diag: dvmul(B->Pmat->info.diag.ndiag,B->Pmat->info.diag.idiag,1,r,1,z,1); break; case Pre_Block:{ register int i; int eDIM = U->fhead->dim(); int nedge = B->Pmat->info.block.nedge; int info,l,cnt; double *tmp = dvector(0,LGmax-1); static double *pinv; dcopy(B->nsolve,r,1,z,1); #ifdef BVERT dpptrs('U',B->Pmat->info.block.nvert,1,B->Pmat->info.block.ivert, z,B->Pmat->info.block.nvert,info); #else dvmul(B->Pmat->info.block.nvert,B->Pmat->info.block.ivert,1,z,1,z,1); #endif cnt = B->Pmat->info.block.nvert; for(i = 0; i < nedge; ++i){ l = B->Pmat->info.block.Ledge[i]; dpptrs('U',l,1,B->Pmat->info.block.iedge[i],z+cnt,l,info); cnt += l; } if(!pinv){ pinv = dvector(0,B->nsolve-cnt-1); for(i = 0; i < B->nsolve-cnt; ++i) pinv[i] = U->flist[i]->get_1diag_massmat(0); dvrecp(B->nsolve-cnt,pinv,1,pinv,1); } /* finally multiply pressure dof by inverse of mass matrix */ dvmul(B->nsolve-cnt,pinv,1,r+cnt,1,z+cnt,1); free(tmp); } break; case Pre_None: dcopy(B->nsolve,r,1,z,1); break; case Pre_Diag_Stokes:{ register int i,j,k; double *tmp = dvector(0,B->nglobal-1),*sc,*b,*p,*wk; int eDIM = U->fhead->dim(); int **bmap = B->bmap; int nel = B->nel,asize,Nbmodes,info; if(eDIM == 2) wk = dvector(0,eDIM*4*LGmax); else wk = dvector(0,eDIM*6*LGmax*LGmax); dzero(B->nglobal,tmp,1); dcopy(B->nsolve-nel,r,1,tmp,1); /* multiply by invc */ dvmul(B->nsolve-nel,B->Pmat->info.diagst.idiag,1,tmp,1,tmp,1); /* multiply by b^T */ sc = B->signchange; p = z + B->nsolve-nel; dcopy(nel,r+B->nsolve-nel,1,p,1); dneg (nel-1,p+1,1); /* start at 1 to omit singular pressure point */ if(!(B->singular)) error_msg(issue in setting singular modes in Stokes_Precon not resolved); for(j = 1; j < nel; ++j){ Nbmodes = U->flist[j]->Nbmodes; asize = Nbmodes*eDIM+1; dzero (asize,wk,1); dgathr(eDIM*Nbmodes, tmp, bmap[j], wk); dvmul (eDIM*Nbmodes, sc, 1, wk, 1, wk, 1); b = B->Gmat->a[U->flist[j]->geom->id] + asize*(asize-1)/2; p[j] += ddot(asize-1,b,1,wk,1); sc += asize; } /* solve for p */ dpptrs('L', nel, 1, B->Pmat->info.diagst.binvcb,p,nel,info); /* generate initial vector as tmp = B p */ sc = B->signchange; dzero(B->nglobal,tmp,1); /* start from 1 due to singular pressure system */ for(k = 1; k < nel; ++k){ Nbmodes = U->flist[k]->Nbmodes; asize = Nbmodes*eDIM+1; b = B->Gmat->a[U->flist[k]->geom->id] + asize*(asize-1)/2; dsmul(asize-1,p[k],b,1,wk,1); for(i = 0; i < eDIM*Nbmodes; ++i) tmp[bmap[k][i]] += sc[i]*wk[i]; sc += asize; } /* set z = r - tmp */ dvsub(B->nsolve-nel,r,1,tmp,1,z,1); /* u = invc*z */ dvmul(B->nsolve-nel,B->Pmat->info.diagst.idiag,1,z,1,z,1); free(tmp); free(wk); } break; case Pre_Block_Stokes:{ register int i,j,k; double *tmp = dvector(0,B->nglobal-1),*sc,*b,*p,*wk; int eDIM = U->fhead->dim(),cnt,l; int **bmap = B->bmap; int nel = B->nel,asize,Nbmodes,info; int nedge = B->Pmat->info.blockst.nedge; if(eDIM == 2) wk = dvector(0,eDIM*4*LGmax); else wk = dvector(0,eDIM*6*LGmax*LGmax); dzero(B->nglobal,tmp,1); dcopy(B->nsolve-nel,r,1,tmp,1); /* multiply by invc */ #ifdef BVERT dpptrs('U',B->Pmat->info.blockst.nvert,1,B->Pmat->info.blockst.ivert, tmp,B->Pmat->info.blockst.nvert,info); #else dvmul(B->Pmat->info.blockst.nvert,B->Pmat->info.blockst.ivert,1, tmp,1,tmp,1); #endif cnt = B->Pmat->info.blockst.nvert; for(i = 0; i < nedge; ++i){ l = B->Pmat->info.blockst.Ledge[i]; dpptrs('U',l,1,B->Pmat->info.blockst.iedge[i],tmp+cnt,l,info); cnt += l; } /* multiply by b^T */ sc = B->signchange; p = z + B->nsolve-nel; dcopy(nel,r+B->nsolve-nel,1,p,1); dneg(nel-1,p+1,1); if(!(B->singular)) error_msg(issue in setting singular modes in Stokes_Precon not resolved); /* start at 1 to omit singular pressure point */ for(j = 1; j < nel; ++j){ Nbmodes = U->flist[j]->Nbmodes; asize = Nbmodes*eDIM+1; dzero (asize,wk,1); dgathr(eDIM*Nbmodes, tmp, bmap[j], wk); dvmul (eDIM*Nbmodes, sc, 1, wk, 1, wk, 1); b = B->Gmat->a[U->flist[j]->geom->id] + asize*(asize-1)/2; p[j] += ddot(asize-1,b,1,wk,1); sc += asize; } /* solve for p */ dpptrs('L', nel, 1, B->Pmat->info.blockst.binvcb,p,nel,info); /* generate initial vector as tmp = B p */ sc = B->signchange; dzero(B->nglobal,tmp,1); /* start from 1 due to singular pressure system */ for(k = 1; k < nel; ++k){ Nbmodes = U->flist[k]->Nbmodes; asize = Nbmodes*eDIM+1; b = B->Gmat->a[U->flist[k]->geom->id] + asize*(asize-1)/2; dsmul(asize-1,p[k],b,1,wk,1); for(i = 0; i < eDIM*Nbmodes; ++i) tmp[bmap[k][i]] += sc[i]*wk[i]; sc += asize; } /* set z = r - tmp */ dvsub(B->nsolve-nel,r,1,tmp,1,z,1); /* u = invc*z */ #ifdef BVERT dpptrs('U',B->Pmat->info.blockst.nvert,1,B->Pmat->info.blockst.ivert, z,B->Pmat->info.blockst.nvert,info); #else dvmul(B->Pmat->info.blockst.nvert,B->Pmat->info.blockst.ivert,1,z,1,z,1); #endif cnt = B->Pmat->info.blockst.nvert; for(i = 0; i < nedge; ++i){ l = B->Pmat->info.blockst.Ledge[i]; dpptrs('U',l,1,B->Pmat->info.blockst.iedge[i],z+cnt,l,info); cnt += l; } free(tmp); free(wk); } break; default: fprintf(stderr,"Preconditioner (Precon=%d) not known\n",B->Precon); exit(-1); } }
static void divergence_free_init(Element_List **V, double *u0, double *r, Bsystem **B, double **wk){ register int i,j,k,cnt; int eDIM = V[0]->fhead->dim(); int asize,info,Nbmodes; int **bmap = B[eDIM]->bmap; int nel = B[eDIM]->nel; int nglobal = B[eDIM]->nglobal; int nsolve = B[eDIM]->nsolve; int vsolve = B[0]->nsolve; double **a = B[eDIM]->Gmat->a,*x,*b,*sc,*sc1; static double *invb; x = dvector(0,nel-1); /* form and invert B' B system */ if(!invb){ invb = dvector(0,nel*(nel+1)/2-1); sc = B[0]->signchange; for(cnt = 0,k = 0; k < nel; ++k){ /* gather B from local systems */ Nbmodes = V[0]->flist[k]->Nbmodes; asize = Nbmodes*eDIM+1; b = a[V[0]->flist[k]->geom->id] + asize*(asize-1)/2; dzero(nglobal,u0,1); for(j = 0; j < eDIM; ++j){ for(i = 0; i < Nbmodes; ++i) u0[bmap[k][i+j*Nbmodes]] += sc[i]*b[i+j*Nbmodes]; } dzero(nglobal-nsolve,u0+nsolve,1); /* take inner product with B' */ sc1 = B[0]->signchange; for(j = k; j < nel; ++j,cnt++){ dzero(asize,wk[0],1); dgathr(asize-1, u0, bmap[j], wk[0]); for(i = 0; i < eDIM; ++i) dvmul (Nbmodes, sc1, 1, wk[0]+i*Nbmodes, 1, wk[0]+i*Nbmodes, 1); Nbmodes = V[0]->flist[j]->Nbmodes; asize = Nbmodes*eDIM+1; b = a[V[0]->flist[j]->geom->id] + asize*(asize-1)/2; invb[cnt] = ddot(asize-1,b,1,wk[0],1); sc1 += Nbmodes; } sc += V[0]->flist[k]->Nbmodes; } /* take out first row to deal with singularity */ if(B[eDIM]->singular) dzero(nel,invb,1); invb[0] = 1.0; dpptrf('L', nel, invb, info); } /* solve B' B x = r */ dcopy (nel,r+eDIM*vsolve,1,x,1); dpptrs('L', nel, 1, invb, x, nel, info); dzero(B[eDIM]->nglobal,u0,1); /* generate initial vector as u0 = B x */ sc = B[0]->signchange; for(k = 0; k < nel; ++k){ Nbmodes = V[0]->flist[k]->Nbmodes; asize = Nbmodes*eDIM+1; b = a[V[0]->flist[k]->geom->id] + asize*(asize-1)/2; dsmul(asize-1,x[k],b,1,wk[0],1); for(j = 0; j < eDIM; ++j){ for(i = 0; i < Nbmodes; ++i) u0[bmap[k][i+j*Nbmodes]] += sc[i]*wk[0][i+j*Nbmodes]; } sc += Nbmodes; } dzero(nglobal-eDIM*vsolve, u0 + eDIM*vsolve, 1); /* subtract off a*u0 from r */ sc = B[0]->signchange; for(k = 0; k < nel; ++k){ Nbmodes = V[0]->flist[k]->Nbmodes; dzero(eDIM*Nbmodes+1,wk[0],1); dgathr(eDIM*Nbmodes+1, u0, bmap[k], wk[0]); for(j = 0; j < eDIM; ++j) dvmul (Nbmodes, sc, 1, wk[0] + j*Nbmodes, 1, wk[0] + j*Nbmodes,1); dspmv('U',eDIM*Nbmodes+1,1.0,a[V[0]->flist[k]->geom->id], wk[0],1,0.0,wk[1],1); for(j = 0; j < eDIM; ++j) for(i = 0; i < Nbmodes; ++i) r[bmap[k][i+j*Nbmodes]] -= sc[i]*wk[1][i+j*Nbmodes]; r[bmap[k][eDIM*Nbmodes]] -= wk[1][eDIM*Nbmodes]; sc += Nbmodes; } r[eDIM*vsolve] = 0.0; dzero(nglobal-nsolve, r + nsolve, 1); free(x); }