static void solve_velocity_interior(Element_List **V, Bsystem **Vbsys){ register int i; int N,nbl,rows,id,info; int *bw = Vbsys[0]->Gmat->bwidth_c; int eDIM = V[0]->flist[0]->dim(); double *hj; int **ipiv = Vbsys[0]->Gmat->cipiv; double **invc = Vbsys[0]->Gmat-> invc; double **binvc = Vbsys[0]->Gmat->binvc; double **invcd = Vbsys[0]->Gmat->invcd; double ***dbinvc = Vbsys[0]->Gmat->dbinvc; Element *E,*P; for(i = 0; i < eDIM; ++i) for(E=V[i]->fhead;E;E=E->next){ N = E->Nmodes - E->Nbmodes; E->state = 't'; if(!N) continue; id = E->geom->id; nbl = E->Nbmodes; hj = E->vert->hj + nbl; P = V[eDIM]->flist[E->id]; rows = P->Nmodes; if(Vbsys[0]->lambda->wave){ /* dbinvc only store dbi in this formulation */ dgemv('N', N, rows, -1., dbinvc[i][id], N, P->vert->hj, 1, 1.0,hj,1); dgetrs('T',N,rows,invc[id],N,ipiv[id],dbinvc[i][id],N,info); if(N > 3*bw[id]) dgbtrs('N', N, bw[id]-1,bw[id]-1, 1, invc[id], 3*bw[id]-2, ipiv[id], hj, N, info); else dgetrs('N', N, 1, invc[id], N, ipiv[id], hj, N, info); dgemv('N', N, nbl , -1., invcd [id], N, E->vert->hj, 1, 1.0,hj,1); } else{ if(N > 2*bw[id]) dpbtrs('L', N, bw[id]-1, 1, invc[id], bw[id], hj, N, info); else dpptrs('L', N, 1, invc[id], hj, N, info); dgemv('N', N, nbl , -1., binvc [id], N, E->vert->hj, 1, 1.0,hj,1); dgemv('N', N, rows, -1., dbinvc[i][id], N, P->vert->hj, 1, 1.0,hj,1); } } }
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); }
static int DTPUMatSolve(void* AA, double b[], double x[], int n){ dtpumat* A=(dtpumat*) AA; ffinteger INFO,NRHS=1,LDB=A->n,N=A->n; double *AP=A->val,*ss=A->sscale; char UPLO=A->UPLO; if (N>0) LDB=N; dtpuscalevec(1.0,ss,b,x,n); dpptrs(&UPLO, &N, &NRHS, AP, x, &LDB, &INFO ); dtpuscalevec(1.0,ss,x,x,n); return INFO; }
void Condense_Recur(Bsystem *Bsys, int lev, char trip){ register int i,j,k; Recur *rdata = Bsys->rslv->rdata + lev; int npatch = rdata->npatch, *bwidth_c = rdata->bwidth_c; int csize,asize,bw,info, *pivot; double *A, *B, *C; double *c, **b; for(i = 0; i < npatch; ++i){ bw = bwidth_c[i] = bandwidth_c(rdata->invc[i],rdata->patchlen_c[i]); asize = rdata->patchlen_a[i]; csize = rdata->patchlen_c[i]; A = Bsys->rslv->A.a[i]; B = rdata->binvc[i]; C = rdata->invc[i]; /* pack and factor C */ if(csize){ if(2*bw < csize) c = dvector(0,csize*bw-1); else{ c = dvector(0,csize*(csize+1)/2-1); if(trip != 'p') pivot = ivector(0,csize-1); } PackMatrixV (C,csize,c,bw,'l'); if(trip == 'p') /* positive definite system */ FacMatrix (c,csize,bw); else{ FacMatrixSym(c,pivot,csize,bw); rdata->pivotc[i] = pivot; } free(rdata->invc[i]); rdata->invc[i] = C = c; if(asize){ /* copy B into b for A-BCB */ b = dmatrix(0,asize-1,0,csize-1); dcopy(asize*csize,B,1,*b,1); /* calc inv(c)*b^T */ if(trip == 'p'){ if(2*bw < csize) dpbtrs('L',csize,bw-1,asize,C, bw,B,csize,info); else dpptrs('L',csize,asize,C,B,csize,info); } else{ if(2*bw < csize){ error_msg(error in codense_recur); } else dsptrs('L',csize,asize,C,pivot,B,csize,info); } /* A - b inv(c) b^T */ for(j = 0; j < asize; ++j) for(k = 0; k < asize; ++k) A[j*asize+k] -= ddot(csize,b[j],1,B+k*csize,1); free_dmatrix(b,0,0); } } } }
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); }