Actual source code: cholbs.c
1: /*$Id: cholbs.c,v 1.64 2001/08/07 03:02:51 balay Exp $*/
3: #include petsc.h
5: /* We must define MLOG for BlockSolve logging */
6: #if defined(PETSC_USE_LOG)
7: #define MLOG
8: #endif
10: #include src/mat/impls/rowbs/mpi/mpirowbs.h
14: int MatCholeskyFactorNumeric_MPIRowbs(Mat mat,Mat *factp)
15: {
16: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)mat->data;
18: #if defined(PETSC_USE_LOG)
19: PetscReal flop1 = BSlocal_flops();
20: #endif
23: /* Do prep work if same nonzero structure as previously factored matrix */
24: if (mbs->factor == FACTOR_CHOLESKY) {
25: /* Copy the nonzeros */
26: BScopy_nz(mbs->pA,mbs->fpA);CHKERRBS(0);
27: }
28: /* Form incomplete Cholesky factor */
29: mbs->0; mbs->failures = 0; mbs->alpha = 1.0;
30: while ((mbs->BSfactor(mbs->fpA,mbs->comm_fpA,mbs->procinfo))) {
31: CHKERRBS(0); mbs->failures++;
32: /* Copy only the nonzeros */
33: BScopy_nz(mbs->pA,mbs->fpA);CHKERRBS(0);
34: /* Increment the diagonal shift */
35: mbs->alpha += 0.1;
36: BSset_diag(mbs->fpA,mbs->alpha,mbs->procinfo);CHKERRBS(0);
37: PetscLogInfo(mat,"MatCholeskyFactorNumeric_MPIRowbs:BlockSolve95: %d failed factor(s), err=%d, alpha=%g\n",
38: mbs->failures,mbs->ierr,mbs->alpha);
39: }
40: #if defined(PETSC_USE_LOG)
41: PetscLogFlops((int)(BSlocal_flops()-flop1));
42: #endif
44: mbs->factor = FACTOR_CHOLESKY;
45: return(0);
46: }
50: int MatLUFactorNumeric_MPIRowbs(Mat mat,Mat *factp)
51: {
52: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)mat->data;
54: #if defined(PETSC_USE_LOG)
55: PetscReal flop1 = BSlocal_flops();
56: #endif
59: /* Do prep work if same nonzero structure as previously factored matrix */
60: if (mbs->factor == FACTOR_LU) {
61: /* Copy the nonzeros */
62: BScopy_nz(mbs->pA,mbs->fpA);CHKERRBS(0);
63: }
64: /* Form incomplete Cholesky factor */
65: mbs->0; mbs->failures = 0; mbs->alpha = 1.0;
66: while ((mbs->BSfactor(mbs->fpA,mbs->comm_fpA,mbs->procinfo))) {
67: CHKERRBS(0); mbs->failures++;
68: /* Copy only the nonzeros */
69: BScopy_nz(mbs->pA,mbs->fpA);CHKERRBS(0);
70: /* Increment the diagonal shift */
71: mbs->alpha += 0.1;
72: BSset_diag(mbs->fpA,mbs->alpha,mbs->procinfo);CHKERRBS(0);
73: PetscLogInfo(mat,"MatLUFactorNumeric_MPIRowbs:BlockSolve95: %d failed factor(s), err=%d, alpha=%g\n",
74: mbs->failures,mbs->ierr,mbs->alpha);
75: }
76: mbs->factor = FACTOR_LU;
77: (*factp)->assembled = PETSC_TRUE;
78: #if defined(PETSC_USE_LOG)
79: PetscLogFlops((int)(BSlocal_flops()-flop1));
80: #endif
81: return(0);
82: }
83: /* ------------------------------------------------------------------- */
86: int MatSolve_MPIRowbs(Mat mat,Vec x,Vec y)
87: {
88: Mat submat = (Mat) mat->data;
89: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)submat->data;
90: int ierr;
91: PetscScalar *ya,*xa,*xworka;
93: #if defined(PETSC_USE_LOG)
94: PetscReal flop1 = BSlocal_flops();
95: #endif
98: /* Permute and apply diagonal scaling to vector, where D^{-1/2} is stored */
99: if (!mbs->vecs_permscale) {
100: VecGetArray(x,&xa);
101: VecGetArray(mbs->xwork,&xworka);
102: BSperm_dvec(xa,xworka,mbs->pA->perm);CHKERRBS(0);
103: VecRestoreArray(x,&xa);
104: VecRestoreArray(mbs->xwork,&xworka);
105: VecPointwiseMult(mbs->diag,mbs->xwork,y);
106: } else {
107: VecCopy(x,y);
108: }
110: VecGetArray(y,&ya);
111: if (mbs->procinfo->single) {
112: /* Use BlockSolve routine for no cliques/inodes */
113: BSfor_solve1(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
114: BSback_solve1(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
115: } else {
116: BSfor_solve(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
117: BSback_solve(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
118: }
119: VecRestoreArray(y,&ya);
121: /* Apply diagonal scaling and unpermute, where D^{-1/2} is stored */
122: if (!mbs->vecs_permscale) {
123: VecPointwiseMult(y,mbs->diag,mbs->xwork);
124: VecGetArray(y,&ya);
125: VecGetArray(mbs->xwork,&xworka);
126: BSiperm_dvec(xworka,ya,mbs->pA->perm);CHKERRBS(0);
127: VecRestoreArray(y,&ya);
128: VecRestoreArray(mbs->xwork,&xworka);
129: }
130: #if defined(PETSC_USE_LOG)
131: PetscLogFlops((int)(BSlocal_flops()-flop1));
132: #endif
133: return(0);
134: }
136: /* ------------------------------------------------------------------- */
139: int MatForwardSolve_MPIRowbs(Mat mat,Vec x,Vec y)
140: {
141: Mat submat = (Mat) mat->data;
142: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)submat->data;
143: int ierr;
144: PetscScalar *ya,*xa,*xworka;
146: #if defined(PETSC_USE_LOG)
147: PetscReal flop1 = BSlocal_flops();
148: #endif
151: /* Permute and apply diagonal scaling to vector, where D^{-1/2} is stored */
152: if (!mbs->vecs_permscale) {
153: VecGetArray(x,&xa);
154: VecGetArray(mbs->xwork,&xworka);
155: BSperm_dvec(xa,xworka,mbs->pA->perm);CHKERRBS(0);
156: VecRestoreArray(x,&xa);
157: VecRestoreArray(mbs->xwork,&xworka);
158: VecPointwiseMult(mbs->diag,mbs->xwork,y);
159: } else {
160: VecCopy(x,y);
161: }
163: VecGetArray(y,&ya);
164: if (mbs->procinfo->single){
165: /* Use BlockSolve routine for no cliques/inodes */
166: BSfor_solve1(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
167: } else {
168: BSfor_solve(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
169: }
170: VecRestoreArray(y,&ya);
172: #if defined(PETSC_USE_LOG)
173: PetscLogFlops((int)(BSlocal_flops()-flop1));
174: #endif
176: return(0);
177: }
179: /* ------------------------------------------------------------------- */
182: int MatBackwardSolve_MPIRowbs(Mat mat,Vec x,Vec y)
183: {
184: Mat submat = (Mat) mat->data;
185: Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)submat->data;
186: int ierr;
187: PetscScalar *ya,*xworka;
189: #if defined (PETSC_USE_LOG)
190: PetscReal flop1 = BSlocal_flops();
191: #endif
194: VecCopy(x,y);
196: VecGetArray(y,&ya);
197: if (mbs->procinfo->single) {
198: /* Use BlockSolve routine for no cliques/inodes */
199: BSback_solve1(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
200: } else {
201: BSback_solve(mbs->fpA,ya,mbs->comm_pA,mbs->procinfo);CHKERRBS(0);
202: }
203: VecRestoreArray(y,&ya);
205: /* Apply diagonal scaling and unpermute, where D^{-1/2} is stored */
206: if (!mbs->vecs_permscale) {
207: VecPointwiseMult(y,mbs->diag,mbs->xwork);
208: VecGetArray(y,&ya);
209: VecGetArray(mbs->xwork,&xworka);
210: BSiperm_dvec(xworka,ya,mbs->pA->perm);CHKERRBS(0);
211: VecRestoreArray(y,&ya);
212: VecRestoreArray(mbs->xwork,&xworka);
213: }
214: #if defined (PETSC_USE_LOG)
215: PetscLogFlops((int)(BSlocal_flops()-flop1));
216: #endif
217: return(0);
218: }
221: /*
222: The logging variables required by BlockSolve,
224: This is an ugly hack that allows PETSc to run properly with BlockSolve regardless
225: of whether PETSc or BlockSolve is compiled with logging turned on.
227: It is bad because it relys on BlockSolve's internals not changing related to
228: logging but we have no choice, plus it is unlikely BlockSolve will be developed
229: in the near future anyways.
230: */
231: double MLOG_flops;
232: double MLOG_event_flops;
233: double MLOG_time_stamp;
234: int MLOG_sequence_num;
235: #if defined (MLOG_MAX_EVNTS)
236: MLOG_log_type MLOG_event_log[MLOG_MAX_EVNTS];
237: MLOG_log_type MLOG_accum_log[MLOG_MAX_ACCUM];
238: #else
239: typedef struct __MLOG_log_type {
240: double time_stamp;
241: double total_time;
242: double flops;
243: int event_num;
244: } MLOG_log_type;
245: #define MLOG_MAX_EVNTS 1300
246: #define MLOG_MAX_ACCUM 75
247: MLOG_log_type MLOG_event_log[MLOG_MAX_EVNTS];
248: MLOG_log_type MLOG_accum_log[MLOG_MAX_ACCUM];
249: #endif