int CleanupSolution(AutoData *Data) { int i; int result = 1; if (Data->u != NULL) { FREE_DMATRIX(Data->u); Data->u = NULL; } if (Data->usm != NULL) { if (Data->usm[0] != NULL) { FREE_DMATRIX(Data->usm[0]); Data->usm[0] = NULL; } if (Data->usm[1] != NULL) { FREE_DMATRIX(Data->usm[1]); Data->usm[1] = NULL; } for (i=0; i<(int)(log2(Data->nsm)); i++) { FREE_DMATRIX(Data->usm[2+i]); Data->usm[2+i] = NULL; } FREE(Data->usm); Data->usm = NULL; } if (Data->par != NULL) { FREE_DMATRIX(Data->par); Data->par = NULL; } if (Data->ev != NULL) { FREE_DCMATRIX(Data->ev); Data->ev = NULL; } if (Data->c0 != NULL) { FREE_DMATRIX_3D(Data->c0); Data->c0 = NULL; } if (Data->c1 != NULL) { FREE_DMATRIX_3D(Data->c1); Data->c1 = NULL; } if (Data->a1 != NULL) { FREE_DMATRIX_3D(Data->a1); Data->a1 = NULL; } if (Data->a2 != NULL) { FREE_DMATRIX_3D(Data->a2); Data->a2 = NULL; } return result; }
int CleanupSpecialPoints(AutoData *Data) { AutoSPData *SPData = Data->sp; int result = 1; int i; if (SPData != NULL) { for (i = 0; i < Data->sp_len; i++) { if (SPData[i].icp != NULL) { FREE(SPData[i].icp); SPData[i].icp = NULL; } if (SPData[i].u != NULL) { FREE(SPData[i].u); SPData[i].u = NULL; } if (SPData[i].rldot != NULL) { FREE(SPData[i].rldot); SPData[i].rldot = NULL; } if (SPData[i].ups != NULL) { FREE_DMATRIX(SPData[i].ups); SPData[i].ups = NULL; } if (SPData[i].udotps != NULL) { FREE_DMATRIX(SPData[i].udotps); SPData[i].udotps = NULL; } if (SPData[i].a1 != NULL) { FREE_DMATRIX_3D(SPData[i].a1); SPData[i].a1 = NULL; } if (SPData[i].a2 != NULL) { FREE_DMATRIX_3D(SPData[i].a2); SPData[i].a2 = NULL; } } FREE(Data->sp); Data->sp = NULL; } return result; }
void FREE_DMATRIX_3D(doublereal ***m) { if (m==NULL) return; FREE_DMATRIX(m[0]); FREE( (char *) (m) ); }
int conpar_mpi_wrapper(integer *nov, integer *na, integer *nra, integer *nca, doublereal ***a, integer *ncb, doublereal ***b, integer *nbc, integer *nrc, doublereal ***c, doublereal *d, integer *irf, integer *icf) { integer loop_start,loop_end; integer loop_start_tmp,loop_end_tmp; int i,comm_size; int *a_counts,*a_displacements; int *b_counts,*b_displacements; int *c_counts,*c_displacements; int *irf_counts,*irf_displacements; int *icf_counts,*icf_displacements; MPI_Comm_size(MPI_COMM_WORLD,&comm_size); a_counts=(int *)MALLOC(sizeof(int)*comm_size); a_displacements=(int *)MALLOC(sizeof(int)*comm_size); b_counts=(int *)MALLOC(sizeof(int)*comm_size); b_displacements=(int *)MALLOC(sizeof(int)*comm_size); c_counts=(int *)MALLOC(sizeof(int)*comm_size); c_displacements=(int *)MALLOC(sizeof(int)*comm_size); irf_counts=(int *)MALLOC(sizeof(int)*comm_size); irf_displacements=(int *)MALLOC(sizeof(int)*comm_size); icf_counts=(int *)MALLOC(sizeof(int)*comm_size); icf_displacements=(int *)MALLOC(sizeof(int)*comm_size); a_counts[0] = 0; a_displacements[0] = 0; b_counts[0] = 0; b_displacements[0] = 0; c_counts[0] = 0; c_displacements[0] = 0; irf_counts[0] = 0; irf_displacements[0] = 0; icf_counts[0] = 0; icf_displacements[0] = 0; for(i=1;i<comm_size;i++){ /*Send message to get worker into conpar mode*/ { int message=AUTO_MPI_CONPAR_MESSAGE; MPI_Send(&message,1,MPI_INT,i,0,MPI_COMM_WORLD); } loop_start = ((i-1)*(*na))/(comm_size - 1); loop_end = ((i)*(*na))/(comm_size - 1); a_counts[i] = (*nca)*(*nra)*(loop_end-loop_start); a_displacements[i] = (*nca)*(*nra)*loop_start; b_counts[i] = (*ncb)*(*nra)*(loop_end-loop_start); b_displacements[i] = (*ncb)*(*nra)*loop_start; c_counts[i] = (*nca)*(*nrc)*(loop_end-loop_start); c_displacements[i] = (*nca)*(*nrc)*loop_start; irf_counts[i] = (*nra)*(loop_end-loop_start); irf_displacements[i] = (*nra)*loop_start; icf_counts[i] = (*nca)*(loop_end-loop_start); icf_displacements[i] = (*nca)*loop_start; loop_start_tmp = 0; loop_end_tmp = loop_end-loop_start; MPI_Send(&loop_start_tmp ,1,MPI_LONG,i,0,MPI_COMM_WORLD); MPI_Send(&loop_end_tmp ,1,MPI_LONG,i,0,MPI_COMM_WORLD); } { integer params[6]; params[0]=*nov; params[1]=*nra; params[2]=*nca; params[3]=*ncb; params[4]=*nbc; params[5]=*nrc; MPI_Bcast(params ,6,MPI_LONG,0,MPI_COMM_WORLD); } MPI_Scatterv(irf,irf_counts,irf_displacements,MPI_LONG, NULL,0,MPI_LONG, 0,MPI_COMM_WORLD); MPI_Scatterv(icf,icf_counts,icf_displacements,MPI_LONG, NULL,0,MPI_LONG, 0,MPI_COMM_WORLD); /* Worker is running now */ { /*I create a temporary send buffer for the MPI_Reduce command. This is because there isn't an asymmetric version (like MPI_Scatterv).*/ double **dtemp; dtemp = DMATRIX(*nrc,*ncb); for(i=0;i<(*nrc);i++) for(j=0;i<(*ncb);i++) dtemp[i][j]=d[i][j]; MPI_Reduce(dtemp[0],d[0],(*ncb)*(*nrc),MPI_DOUBLE,MPI_SUM,0,MPI_COMM_WORLD); FREE_DMATRIX(dtemp); } MPI_Gatherv(NULL,0,MPI_DOUBLE, a[0][0],a_counts,a_displacements,MPI_DOUBLE, 0,MPI_COMM_WORLD); MPI_Gatherv(NULL,0,MPI_DOUBLE, b[0][0],b_counts,b_displacements,MPI_DOUBLE, 0,MPI_COMM_WORLD); MPI_Gatherv(NULL,0,MPI_DOUBLE, c[0][0],c_counts,c_displacements,MPI_DOUBLE, 0,MPI_COMM_WORLD); MPI_Gatherv(NULL,0,MPI_LONG, irf,irf_counts,irf_displacements,MPI_LONG, 0,MPI_COMM_WORLD); MPI_Gatherv(NULL,0,MPI_LONG, icf,icf_counts,icf_displacements,MPI_LONG, 0,MPI_COMM_WORLD); return 0; }
void *conpar_process(void * arg) { integer icf_dim1, irf_dim1; /* Local variables */ integer ipiv, jpiv, itmp; doublereal tpiv; integer i, j, l, k1, k2, m1, m2, ic, ir; doublereal rm; integer ir1, irp; doublereal piv; integer icp1; integer *nov, *nra, *nca; doublereal ***a; integer *ncb; doublereal ***b; integer *nbc, *nrc; doublereal ***c, **d; integer *irf, *icf; integer loop_start,loop_end; #ifdef PTHREADS doublereal **dlocal; #endif #ifdef USAGE struct rusage *conpar_process_usage; usage_start(&conpar_process_usage); #endif nov = ((conpar_parallel_arglist *)arg)->nov; nra = ((conpar_parallel_arglist *)arg)->nra; nca = ((conpar_parallel_arglist *)arg)->nca; a = ((conpar_parallel_arglist *)arg)->a; ncb = ((conpar_parallel_arglist *)arg)->ncb; b = ((conpar_parallel_arglist *)arg)->b; nbc = ((conpar_parallel_arglist *)arg)->nbc; nrc = ((conpar_parallel_arglist *)arg)->nrc; c = ((conpar_parallel_arglist *)arg)->c; d = ((conpar_parallel_arglist *)arg)->d; irf = ((conpar_parallel_arglist *)arg)->irf; icf = ((conpar_parallel_arglist *)arg)->icf; loop_start = ((conpar_parallel_arglist *)arg)->loop_start; loop_end = ((conpar_parallel_arglist *)arg)->loop_end; #ifdef PTHREADS dlocal=DMATRIX(*nrc, *ncb); #endif /* In the default case we don't need to do anything special */ if(global_conpar_type == CONPAR_DEFAULT) { ; } /* In the message passing case we set d to be 0.0, do a sum here, and then do the final sum (with the true copy of d) in the master */ else if (global_conpar_type == CONPAR_MPI) { for(i=0;i<*nrc;i++) for (j=0; j<*ncb;j++) d[i][j]=0.0; } /* In the shared memory case we create a local variable for doing this threads part of the sum, then we do a final sum into shared memory at the end */ else if (global_conpar_type == CONPAR_PTHREADS) { #ifdef PTHREADS for(i=0;i<*nrc;i++) for (j=0; j<*ncb;j++) dlocal[i][j]=0.0; #else ; #endif } /* Note that the summation of the adjacent overlapped part of C */ /* is delayed until REDUCE, in order to merge it with other communications.*/ /* NA is the local NTST. */ irf_dim1 = *nra; icf_dim1 = *nca; /* Condensation of parameters (Elimination of local variables). */ m1 = *nov + 1; m2 = *nca - *nov; for (i = loop_start;i < loop_end; i++) { for (ic = m1; ic <= m2; ++ic) { ir1 = ic - *nov + 1; irp = ir1 - 1; icp1 = ic + 1; /* **Search for pivot (Complete pivoting) */ piv = 0.0; ipiv = irp; jpiv = ic; for (k1 = irp; k1 <= *nra; ++k1) { int irf_k1_i = irf[-1 + k1 + i*irf_dim1]; for (k2 = ic; k2 <= m2; ++k2) { int icf_k2_i = icf[-1 + k2 + i*icf_dim1]; tpiv = a[i][-1 + irf_k1_i][-1 + icf_k2_i]; if (tpiv < 0.0) { tpiv = -tpiv; } if (piv < tpiv) { piv = tpiv; ipiv = k1; jpiv = k2; } } } /* **Move indices */ itmp = icf[-1 + ic + i*icf_dim1]; icf[-1 + ic + i*icf_dim1] = icf[-1 + jpiv + i*icf_dim1]; icf[-1 + jpiv + i*icf_dim1] = itmp; itmp = irf[-1 + irp + i*irf_dim1]; irf[-1 + irp + i*irf_dim1] = irf[-1 + ipiv + i*irf_dim1]; irf[-1 + ipiv + i*irf_dim1] = itmp; { int icf_ic_i = icf[-1 + ic + i*icf_dim1]; int irf_irp_i = irf[-1 + irp + i*irf_dim1]; doublereal *a_offset2 = a[i][-1 + irf_irp_i]; doublereal *b_offset2 = b[i][-1 + irf_irp_i]; /* **End of pivoting; elimination starts here */ for (ir = ir1; ir <= *nra; ++ir) { int irf_ir_i = irf[-1 + ir + i*irf_dim1]; doublereal *a_offset1 = a[i][-1 + irf_ir_i]; doublereal *b_offset1 = b[i][-1 + irf_ir_i]; rm = a_offset1[-1 + icf_ic_i]/a_offset2[-1 + icf_ic_i]; a_offset1[-1 + icf_ic_i] = rm; if (rm != (double)0.) { for (l = 0; l < *nov; ++l) { a_offset1[l] -= rm * a_offset2[l]; } for (l = icp1 -1; l < *nca; ++l) { int icf_l_i = icf[l + i*icf_dim1]; a_offset1[-1 + icf_l_i] -= rm * a_offset2[-1 + icf_l_i]; } for (l = 0; l < *ncb; ++l) { b_offset1[l] -= rm * b_offset2[l]; } } } for (ir = *nbc + 1; ir <= *nrc; ++ir) { doublereal *c_offset1 = c[i][-1 + ir]; doublereal *d_offset1 = d[-1 + ir]; rm = c_offset1[-1 + icf_ic_i]/a_offset2[-1 + icf_ic_i]; c_offset1[-1 + icf_ic_i]=rm; if (rm != (double)0.) { for (l = 0; l < *nov; ++l) { c_offset1[l] -= rm * a_offset2[l]; } for (l = icp1 -1 ; l < *nca; ++l) { int icf_l_i = icf[l + i*icf_dim1]; c_offset1[-1 + icf_l_i] -= rm * a_offset2[-1 + icf_l_i]; } for (l = 0; l < *ncb; ++l) { /* A little explanation of what is going on here is in order I believe. This array is created by a summation across all workers, hence it needs a mutex to avoid concurrent writes (in the shared memory case) or a summation in the master (in the message passing case). Since mutex's can be somewhat slow, we will do the summation into a local variable, and then do a final summation back into global memory when the main loop is done. */ /* Nothing special for the default case */ if(global_conpar_type == CONPAR_DEFAULT) { d_offset1[l] -= rm * b_offset2[l]; } /* In the message passing case we sum into d, which is a local variable initialized to 0.0. We then sum our part with the masters part in the master. */ else if (global_conpar_type == CONPAR_MPI) { d_offset1[l] -= rm * b_offset2[l]; } /* In the shared memory case we sum into a local variable our contribution, and then sum into shared memory at the end (inside a mutex */ else if (global_conpar_type == CONPAR_PTHREADS) { #ifdef PTHREADS dlocal[-1 + ir][l] -= rm * b_offset2[l]; #else ; #endif } } } } } } } #ifdef PTHREADS /* This is were we sum into the global copy of the d array, in the shared memory case */ if(global_conpar_type == CONPAR_PTHREADS) { #ifdef PTHREADS pthread_mutex_lock(&mutex_for_d); for(i=0;i<*nrc;i++) for (j=0; j<*ncb;j++) d[i][j] += dlocal[i][j]; pthread_mutex_unlock(&mutex_for_d); FREE_DMATRIX(dlocal); #else ; #endif } #endif #ifdef USAGE usage_end(conpar_process_usage,"in conpar worker"); #endif return NULL; }
// Assumes there has been a call to SetData before this. int SetInitPoint(double *u, int npar, int *ipar, double *par, int *icp, int nups, double *ups, double *udotps, double *rldot, int adaptcycle) { /* Still think I should get rid of typemap(freearg) and send address to pointer in prepare_cycle so I can realloc if necessary and not have to use this my_... crap */ integer itp; integer i, j; integer SUCCESS = 1; double *my_ups = ups; double *my_udotps = udotps; double *my_rldot = rldot; if (ups == NULL) itp = 3; else { itp = 9; if (adaptcycle || udotps == NULL || rldot == NULL) { // NOTE: AUTO allocates (ntst+1)*ncol points, but only displays ntpl=ntst*ncol+1 my_ups = (double *)MALLOC(gIData->iap.ncol*(gIData->iap.ntst+1)*(gIData->iap.ndim+1)*sizeof(double)); my_udotps = (double *)MALLOC(gIData->iap.ncol*(gIData->iap.ntst+1)*gIData->iap.ndim*sizeof(double)); my_rldot = (double *)MALLOC(gIData->iap.nicp*sizeof(double)); prepare_cycle(gIData, ups, nups, my_ups, my_udotps, my_rldot); } } if (!CreateSpecialPoint(gIData, itp, 1, u, npar, ipar, par, icp, my_ups, my_udotps, my_rldot)) { fprintf(stderr,"*** Warning [interface.c]: Problem in CreateSpecialPoint.\n"); fflush(stderr); SUCCESS = 0; } if ((SUCCESS) && (itp == 9) && (adaptcycle || udotps == NULL || rldot == NULL)) { integer nmx, npr, verbosity; double ds, dsmin, dsmax; // Adjust rldot and sp.nfpr = 1 (AUTO detects this to get udotps and rldot for // starting point) gIData->sp[0].nfpr = 1; // Remove from here till } once findPeriodicOrbit is created FREE(my_ups); FREE(my_udotps); FREE(my_rldot); // Make sure initial point is on curve... nmx = gIData->iap.nmx; npr = gIData->iap.npr; ds = gIData->rap.ds; dsmin = gIData->rap.dsmin; dsmax = gIData->rap.dsmax; verbosity = gIData->verbosity; gIData->iap.nmx = 3; gIData->iap.npr = 3; gIData->rap.ds = min(1e-4, gIData->rap.ds); gIData->rap.dsmin = min(1e-4, gIData->rap.dsmin); gIData->rap.dsmax = min(1e-4, gIData->rap.dsmax); gIData->verbosity = 0; AUTO(gIData); CleanupSolution(gIData); gIData->iap.nmx = nmx; gIData->iap.npr = npr; gIData->rap.ds = ds; gIData->rap.dsmin = dsmin; gIData->rap.dsmax = dsmax; gIData->verbosity = verbosity; // Check for NaNs for (i=0; i<gIData->sp[0].nar; i++) { if (isnan(gIData->sp[1].u[i])) { fprintf(stderr,"*** Warning [interface.c]: NaNs in auto solution.\n"); fflush(stderr); SUCCESS = 0; break; } } if (SUCCESS) { for (i=0; i<gIData->sp[0].nar; i++) gIData->sp[0].u[i] = gIData->sp[1].u[i]; for (i=0; i<gIData->iap.ntst*gIData->iap.ncol+1; i++) for (j=0; j<gIData->iap.ndim+1; j++) gIData->sp[0].ups[i][j] = gIData->sp[1].ups[i][j]; for (i=0; i<gIData->iap.ntst*gIData->iap.ncol+1; i++) for (j=0; j<gIData->iap.ndim; j++) gIData->sp[0].udotps[i][j] = gIData->sp[1].udotps[i][j]; for (i=0; i<gIData->iap.nicp; i++) gIData->sp[0].rldot[i] = gIData->sp[1].rldot[i]; } gIData->sp[0].nfpr = gIData->iap.nicp; FREE(gIData->sp[1].icp); FREE(gIData->sp[1].u); FREE(gIData->sp[1].rldot); FREE_DMATRIX(gIData->sp[1].ups); FREE_DMATRIX(gIData->sp[1].udotps); if (gIData->sflow) { FREE_DMATRIX_3D(gIData->sp[1].a1); FREE_DMATRIX_3D(gIData->sp[1].a2); } gIData->sp[1].icp = NULL; gIData->sp[1].u = NULL; gIData->sp[1].rldot = NULL; gIData->sp[1].ups = NULL; gIData->sp[1].udotps = NULL; gIData->sp[1].a1 = NULL; gIData->sp[1].a2 = NULL; gIData->num_sp = 1; gIData->sp_len = 1; gIData->sp = (AutoSPData *)REALLOC(gIData->sp, (gIData->num_sp)*sizeof(AutoSPData)); } return SUCCESS; }