void FCV_BBDREINIT(long int *Nloc, long int *mudq, long int *mldq, realtype* dqrely, int *ier) { /* First call CVReInitBBD to re-initialize CVBBDPRE module: CVBBD_Data is the pointer to P_data mudq,mldq are the half-bandwidths for computing preconditioner blocks dqrely is the difference quotient relative increment factor FCVgloc is a pointer to the CVLocalFn function FCVcfn is a pointer to the CVCommFn function */ *ier = CVBBDPrecReInit(CVBBD_Data, *mudq, *mldq, *dqrely, FCVgloc, FCVcfn); }
int CVBBDPrecReInitB(void *cvode_mem, int which, int mudqB, int mldqB, realtype dqrelyB) { CVodeMem cv_mem; CVadjMem ca_mem; CVodeBMem cvB_mem; void *cvodeB_mem; int flag; /* Check if cvode_mem exists */ if (cvode_mem == NULL) { cvProcessError(NULL, CVSPILS_MEM_NULL, "CVBBDPRE", "CVBBDPrecReInitB", MSGBBD_MEM_NULL); return(CVSPILS_MEM_NULL); } cv_mem = (CVodeMem) cvode_mem; /* Was ASA initialized? */ if (cv_mem->cv_adjMallocDone == FALSE) { cvProcessError(cv_mem, CVSPILS_NO_ADJ, "CVBBDPRE", "CVBBDPrecReInitB", MSGBBD_NO_ADJ); return(CVSPILS_NO_ADJ); } ca_mem = cv_mem->cv_adj_mem; /* Check which */ if ( which >= ca_mem->ca_nbckpbs ) { cvProcessError(cv_mem, CVSPILS_ILL_INPUT, "CVBBDPRE", "CVBBDPrecReInitB", MSGBBD_BAD_WHICH); return(CVSPILS_ILL_INPUT); } /* Find the CVodeBMem entry in the linked list corresponding to which */ cvB_mem = ca_mem->cvB_mem; while (cvB_mem != NULL) { if ( which == cvB_mem->cv_index ) break; cvB_mem = cvB_mem->cv_next; } cvodeB_mem = (void *) (cvB_mem->cv_mem); flag = CVBBDPrecReInit(cvodeB_mem, mudqB, mldqB, dqrelyB); return(flag); }
int CVBBDPrecReInitB(void *cvadj_mem, long int mudqB, long int mldqB, realtype dqrelyB, CVLocalFnB glocB, CVCommFnB cfnB) { CVadjMem ca_mem; void *cvode_mem; int flag; if (cvadj_mem == NULL) return(CV_ADJMEM_NULL); ca_mem = (CVadjMem) cvadj_mem; cvode_mem = (void *) ca_mem->cvb_mem; gloc_B = glocB; cfn_B = cfnB; flag = CVBBDPrecReInit(bbd_data_B, mudqB, mldqB, dqrelyB, CVAgloc, CVAcfn); return(flag); }
int main(int argc, char *argv[]) { UserData data; void *cvode_mem; realtype abstol, reltol, t, tout; N_Vector u; int iout, my_pe, npes, flag, jpre; long int neq, local_N, mudq, mldq, mukeep, mlkeep; MPI_Comm comm; data = NULL; cvode_mem = NULL; u = NULL; /* Set problem size neq */ neq = NVARS*MX*MY; /* Get processor number and total number of pe's */ MPI_Init(&argc, &argv); comm = MPI_COMM_WORLD; MPI_Comm_size(comm, &npes); MPI_Comm_rank(comm, &my_pe); if (npes != NPEX*NPEY) { if (my_pe == 0) fprintf(stderr, "\nMPI_ERROR(0): npes = %d is not equal to NPEX*NPEY = %d\n\n", npes, NPEX*NPEY); MPI_Finalize(); return(1); } /* Set local length */ local_N = NVARS*MXSUB*MYSUB; /* Allocate and load user data block */ data = (UserData) malloc(sizeof *data); if(check_flag((void *)data, "malloc", 2, my_pe)) MPI_Abort(comm, 1); InitUserData(my_pe, local_N, comm, data); /* Allocate and initialize u, and set tolerances */ u = N_VNew_Parallel(comm, local_N, neq); if(check_flag((void *)u, "N_VNew_Parallel", 0, my_pe)) MPI_Abort(comm, 1); SetInitialProfiles(u, data); abstol = ATOL; reltol = RTOL; /* Call CVodeCreate to create the solver memory and specify the * Backward Differentiation Formula and the use of a Newton iteration */ cvode_mem = CVodeCreate(CV_BDF, CV_NEWTON); if(check_flag((void *)cvode_mem, "CVodeCreate", 0, my_pe)) MPI_Abort(comm, 1); /* Set the pointer to user-defined data */ flag = CVodeSetUserData(cvode_mem, data); if(check_flag(&flag, "CVodeSetUserData", 1, my_pe)) MPI_Abort(comm, 1); /* Call CVodeInit to initialize the integrator memory and specify the * user's right hand side function in u'=f(t,u), the inital time T0, and * the initial dependent variable vector u. */ flag = CVodeInit(cvode_mem, f, T0, u); if(check_flag(&flag, "CVodeInit", 1, my_pe)) return(1); /* Call CVodeSStolerances to specify the scalar relative tolerance * and scalar absolute tolerances */ flag = CVodeSStolerances(cvode_mem, reltol, abstol); if (check_flag(&flag, "CVodeSStolerances", 1, my_pe)) return(1); /* Call CVSpgmr to specify the linear solver CVSPGMR with left preconditioning and the default maximum Krylov dimension maxl */ flag = CVSpgmr(cvode_mem, PREC_LEFT, 0); if(check_flag(&flag, "CVBBDSpgmr", 1, my_pe)) MPI_Abort(comm, 1); /* Initialize BBD preconditioner */ mudq = mldq = NVARS*MXSUB; mukeep = mlkeep = NVARS; flag = CVBBDPrecInit(cvode_mem, local_N, mudq, mldq, mukeep, mlkeep, ZERO, flocal, NULL); if(check_flag(&flag, "CVBBDPrecAlloc", 1, my_pe)) MPI_Abort(comm, 1); /* Print heading */ if (my_pe == 0) PrintIntro(npes, mudq, mldq, mukeep, mlkeep); /* Loop over jpre (= PREC_LEFT, PREC_RIGHT), and solve the problem */ for (jpre = PREC_LEFT; jpre <= PREC_RIGHT; jpre++) { /* On second run, re-initialize u, the integrator, CVBBDPRE, and CVSPGMR */ if (jpre == PREC_RIGHT) { SetInitialProfiles(u, data); flag = CVodeReInit(cvode_mem, T0, u); if(check_flag(&flag, "CVodeReInit", 1, my_pe)) MPI_Abort(comm, 1); flag = CVBBDPrecReInit(cvode_mem, mudq, mldq, ZERO); if(check_flag(&flag, "CVBBDPrecReInit", 1, my_pe)) MPI_Abort(comm, 1); flag = CVSpilsSetPrecType(cvode_mem, PREC_RIGHT); check_flag(&flag, "CVSpilsSetPrecType", 1, my_pe); if (my_pe == 0) { printf("\n\n-------------------------------------------------------"); printf("------------\n"); } } if (my_pe == 0) { printf("\n\nPreconditioner type is: jpre = %s\n\n", (jpre == PREC_LEFT) ? "PREC_LEFT" : "PREC_RIGHT"); } /* In loop over output points, call CVode, print results, test for error */ for (iout = 1, tout = TWOHR; iout <= NOUT; iout++, tout += TWOHR) { flag = CVode(cvode_mem, tout, u, &t, CV_NORMAL); if(check_flag(&flag, "CVode", 1, my_pe)) break; PrintOutput(cvode_mem, my_pe, comm, u, t); } /* Print final statistics */ if (my_pe == 0) PrintFinalStats(cvode_mem); } /* End of jpre loop */ /* Free memory */ N_VDestroy_Parallel(u); free(data); CVodeFree(&cvode_mem); MPI_Finalize(); return(0); }