int new_linearisation(struct jacobian *jac,double hinvGak,int neq,ErrorMsg error_message){ double luparity, *Ax; int i,j,*Ap,*Ai,funcreturn; if(jac->use_sparse==1){ Ap = jac->spJ->Ap; Ai = jac->spJ->Ai; Ax = jac->spJ->Ax; /* Construct jac->spJ->Ax from jac->xjac, the jacobian:*/ for(j=0;j<neq;j++){ for(i=Ap[j];i<Ap[j+1];i++){ if(Ai[i]==j){ /* I'm at the diagonal */ Ax[i] = 1.0-hinvGak*jac->xjac[i]; } else{ Ax[i] = -hinvGak*jac->xjac[i]; } } } /* Matrix constructed... */ if(jac->new_jacobian==_TRUE_){ /*I have a new pattern, and I have not done a LU decomposition since the last jacobian calculation, so I need to do a full sparse LU-decomposition: */ /* Find the sparsity pattern C = J + J':*/ calc_C(jac); /* Calculate the optimal ordering: */ sp_amd(jac->Cp, jac->Ci, neq, jac->cnzmax, jac->Numerical->q,jac->Numerical->wamd); /* if the next line is uncomented, the code uses natural ordering instead of AMD ordering */ /*jac->Numerical->q = NULL;*/ funcreturn = sp_ludcmp(jac->Numerical, jac->spJ, 1e-3); class_test(funcreturn == _FAILURE_,error_message, "Failure in sp_ludcmp. Possibly singular matrix!"); jac->new_jacobian = _FALSE_; } else{ /* I have a repeated pattern, so I can just refactor:*/ sp_refactor(jac->Numerical, jac->spJ); } } else{ /* Normal calculation: */ for(i=1;i<=neq;i++){ for(j=1;j<=neq;j++){ jac->LU[i][j] = - hinvGak * jac->dfdy[i][j]; if(i==j) jac->LU[i][j] +=1.0; } } /*Dense LU decomposition: */ funcreturn = ludcmp(jac->LU,neq,jac->luidx,&luparity,jac->LUw); class_test(funcreturn == _FAILURE_,error_message, "Failure in ludcmp. Possibly singular matrix!"); } return _SUCCESS_; }
int parser_read_double( struct file_content * pfc, char * name, double * value, int * found, ErrorMsg errmsg ) { int index; int i; /* intialize the 'found' flag to false */ * found = _FALSE_; /* search parameter */ index=0; while ((index < pfc->size) && (strcmp(pfc->name[index],name) != 0)) index++; /* if parameter not found, return with 'found' flag still equal to false */ if (index == pfc->size) return _SUCCESS_; /* read parameter value. If this fails, return an error */ class_test(sscanf(pfc->value[index],"%lg",value) != 1, errmsg, "could not read value of parameter %s in file %s\n",name,pfc->filename); /* if parameter read correctly, set 'found' flag to true, as well as the flag associated with this parameter in the file_content structure */ * found = _TRUE_; pfc->read[index] = _TRUE_; /* check for multiple entries of the same parameter. If another occurence is found, return an error. */ for (i=index+1; i < pfc->size; i++) { class_test(strcmp(pfc->name[i],name) == 0, errmsg, "multiple entry of parameter %s in file %s\n",name,pfc->filename); } /* if everything proceeded normally, return with 'found' flag equal to true */ return _SUCCESS_; }
int parser_cat( struct file_content * pfc1, struct file_content * pfc2, struct file_content * pfc3, ErrorMsg errmsg ) { int i; class_test(pfc1->size < 0., errmsg, "size of file_content structure probably not initialized properly\n"); class_test(pfc2->size < 0., errmsg, "size of file_content structure probably not initialized properly\n"); if (pfc1->size == 0) { class_alloc(pfc3->filename,(strlen(pfc2->filename)+1)*sizeof(char),errmsg); sprintf(pfc3->filename,"%s",pfc2->filename); } if (pfc2->size == 0) { class_alloc(pfc3->filename,(strlen(pfc1->filename)+1)*sizeof(char),errmsg); sprintf(pfc3->filename,"%s",pfc1->filename); } if ((pfc1->size !=0) && (pfc2->size != 0)) { class_alloc(pfc3->filename,(strlen(pfc1->filename)+strlen(pfc2->filename)+5)*sizeof(char),errmsg); sprintf(pfc3->filename,"%s or %s",pfc1->filename,pfc2->filename); } pfc3->size = pfc1->size + pfc2->size; class_alloc(pfc3->value,pfc3->size*sizeof(FileArg),errmsg); class_alloc(pfc3->name,pfc3->size*sizeof(FileArg),errmsg); class_alloc(pfc3->read,pfc3->size*sizeof(short),errmsg); for (i=0; i < pfc1->size; i++) { strcpy(pfc3->value[i],pfc1->value[i]); strcpy(pfc3->name[i],pfc1->name[i]); pfc3->read[i]=pfc1->read[i]; } for (i=0; i < pfc2->size; i++) { strcpy(pfc3->value[i+pfc1->size],pfc2->value[i]); strcpy(pfc3->name[i+pfc1->size],pfc2->name[i]); pfc3->read[i+pfc1->size]=pfc2->read[i]; } return _SUCCESS_; }
int parser_read_file( char * filename, struct file_content * pfc, ErrorMsg errmsg ){ FILE * inputfile; char line[_LINE_LENGTH_MAX_]; int counter; int is_data; FileArg name; FileArg value; class_open(inputfile,filename,"r",errmsg); counter = 0; while (fgets(line,_LINE_LENGTH_MAX_,inputfile) != NULL) { class_call(parser_read_line(line,&is_data,name,value,errmsg),errmsg,errmsg); if (is_data == _TRUE_) counter++; } class_test(counter == 0, errmsg, "No readable input in file %s",filename); class_alloc(pfc->filename,(strlen(filename)+1)*sizeof(char),errmsg); strcpy(pfc->filename,filename); class_call(parser_init(pfc,counter,errmsg), errmsg, errmsg); rewind(inputfile); counter = 0; while (fgets(line,_LINE_LENGTH_MAX_,inputfile) != NULL) { class_call(parser_read_line(line,&is_data,name,value,errmsg),errmsg,errmsg); if (is_data == _TRUE_) { strcpy(pfc->name[counter],name); strcpy(pfc->value[counter],value); pfc->read[counter]=_FALSE_; counter++; } } fclose(inputfile); return _SUCCESS_; }
int parser_read_line( char * line, int * is_data, char * name, char * value, ErrorMsg errmsg ) { char * phash; char * pequal; char * left; char * right; /* check that there is an '=' */ pequal=strchr(line,'='); if (pequal == NULL) {*is_data = _FALSE_; return _SUCCESS_;} /* if yes, check that there is not an '#' before the '=' */ phash=strchr(line,'#'); if ((phash != NULL) && (phash-pequal<2)) {*is_data = _FALSE_; return _SUCCESS_;} /* get the name, i.e. the block before the '=' */ left=line; while (left[0]==' ') { left++; } right=pequal-1; while (right[0]==' ') { right--; } if (right-left < 0) {*is_data = _FALSE_; return _SUCCESS_;} class_test(right-left+1 >= _ARGUMENT_LENGTH_MAX_, errmsg, "name starting by '%s' too long; shorten it or increase _ARGUMENT_LENGTH_MAX_",strncpy(name,left,(_ARGUMENT_LENGTH_MAX_-1))); strncpy(name,left,right-left+1); name[right-left+1]='\0'; /* get the value, i.e. the block after the '=' */ left = pequal+1; while (left[0]==' ') { left++; } if (phash == NULL) right = line+strlen(line)-1; else right = phash-1; while (right[0]<=' ') { right--; } if (right-left < 0) {*is_data = _FALSE_; return _SUCCESS_;} class_test(right-left+1 >= _ARGUMENT_LENGTH_MAX_, errmsg, "value starting by '%s' too long; shorten it or increase _ARGUMENT_LENGTH_MAX_",strncpy(value,left,(_ARGUMENT_LENGTH_MAX_-1))); strncpy(value,left,right-left+1); value[right-left+1]='\0'; *is_data = _TRUE_; return _SUCCESS_; }
int parser_read_list_of_strings( struct file_content * pfc, char * name, int * size, char ** pointer_to_list, int * found, ErrorMsg errmsg ) { int index; int i; char * string; char * substring; FileArg string_with_one_value; char * list; /* intialize the 'found' flag to false */ * found = _FALSE_; /* search parameter */ index=0; while ((index < pfc->size) && (strcmp(pfc->name[index],name) != 0)) index++; /* if parameter not found, return with 'found' flag still equal to false */ if (index == pfc->size) return _SUCCESS_; /* count number of comas and compute size = number of comas + 1 */ i = 0; string = pfc->value[index]; do { i ++; substring = strchr(string,','); string = substring+1; } while(substring != NULL); *size = i; /* free and re-allocate array of values */ class_alloc(list,*size*sizeof(FileArg),errmsg); *pointer_to_list = list; /* read one string between each comas */ i = 0; string = pfc->value[index]; do { i ++; substring = strchr(string,','); if (substring == NULL) { strcpy(string_with_one_value,string); } else { strncpy(string_with_one_value,string,(substring-string)); string_with_one_value[substring-string]='\0'; } strcpy(list+(i-1)*_ARGUMENT_LENGTH_MAX_,string_with_one_value); //Insert EOL character: *(list+i*_ARGUMENT_LENGTH_MAX_-1) = '\n'; string = substring+1; } while(substring != NULL); /* if parameter read correctly, set 'found' flag to true, as well as the flag associated with this parameter in the file_content structure */ * found = _TRUE_; pfc->read[index] = _TRUE_; /* check for multiple entries of the same parameter. If another occurence is found, return an error. */ for (i=index+1; i < pfc->size; i++) { class_test(strcmp(pfc->name[i],name) == 0, errmsg, "multiple entry of parameter %s in file %s\n",name,pfc->filename); } /* if everything proceeded normally, return with 'found' flag equal to true */ return _SUCCESS_; }
int evolver_ndf15( int (*derivs)(double x,double * y,double * dy, void * parameters_and_workspace, ErrorMsg error_message), double x_ini, double x_final, double * y_inout, int * used_in_output, int neq, void * parameters_and_workspace_for_derivs, double rtol, double minimum_variation, int (*timescale_and_approximation)(double x, void * parameters_and_workspace, double * timescales, ErrorMsg error_message), double timestep_over_timescale, double * t_vec, int tres, int (*output)(double x,double y[],double dy[],int index_x,void * parameters_and_workspace, ErrorMsg error_message), int (*print_variables)(double x, double y[], double dy[], void *parameters_and_workspace, ErrorMsg error_message), ErrorMsg error_message){ /* Constants: */ double G[5]={1.0,3.0/2.0,11.0/6.0,25.0/12.0,137.0/60.0}; double alpha[5]={-37.0/200,-1.0/9.0,-8.23e-2,-4.15e-2, 0}; double invGa[5],erconst[5]; double abstol = 1e-18, eps=1e-19, threshold=abstol; int maxit=4, maxk=5; /* Logicals: */ int Jcurrent,havrate,done,at_hmin,nofailed,gotynew,tooslow,*interpidx; /* Storage: */ double *f0,*y,*wt,*ddfddt,*pred,*ynew,*invwt,*rhs,*psi,*difkp1,*del,*yinterp; double *tempvec1,*tempvec2,*ypinterp,*yppinterp; double **dif; struct jacobian jac; struct numjac_workspace nj_ws; /* Method variables: */ double t,t0,tfinal,tnew=0; double rh,htspan,absh,hmin,hmax,h,tdel; double abshlast,hinvGak,minnrm,oldnrm=0.,newnrm; double err,hopt,errkm1,hkm1,errit,rate=0.,temp,errkp1,hkp1,maxtmp; int k,klast,nconhk,iter,next,kopt,tdir; /* Misc: */ int stepstat[6],nfenj,j,ii,jj, numidx, neqp=neq+1; int verbose=0; /** Allocate memory . */ void * buffer; int buffer_size; buffer_size= 15*neqp*sizeof(double) +neqp*sizeof(int) +neqp*sizeof(double*) +(7*neq+1)*sizeof(double); class_alloc(buffer, buffer_size, error_message); f0 =(double*)buffer; wt =f0+neqp; ddfddt =wt+neqp; pred =ddfddt+neqp; y =pred+neqp; invwt =y+neqp; rhs =invwt+neqp; psi =rhs+neqp; difkp1 =psi+neqp; del =difkp1+neqp; yinterp =del+neqp; ypinterp =yinterp+neqp; yppinterp=ypinterp+neqp; tempvec1 =yppinterp+neqp; tempvec2 =tempvec1+neqp; interpidx=(int*)(tempvec2+neqp); dif =(double**)(interpidx+neqp); dif[1] =(double*)(dif+neqp); for(j=2;j<=neq;j++) dif[j] = dif[j-1]+7; /* Set row pointers... */ dif[0] = NULL; /* for (ii=0;ii<(7*neq+1);ii++) dif[1][ii]=0.; */ for (j=1; j<=neq; j++) { for (ii=1;ii<=7;ii++) { dif[j][ii]=0.; } } /* class_alloc(f0,sizeof(double)*neqp,error_message); */ /* class_alloc(wt,sizeof(double)*neqp,error_message); */ /* class_alloc(ddfddt,sizeof(double)*neqp,error_message); */ /* class_alloc(pred,sizeof(double)*neqp,error_message); */ /* class_alloc(y,sizeof(double)*neqp,error_message); */ /* class_alloc(invwt,sizeof(double)*neqp,error_message); */ /* class_alloc(rhs,sizeof(double)*neqp,error_message); */ /* class_alloc(psi,sizeof(double)*neqp,error_message); */ /* class_alloc(difkp1,sizeof(double)*neqp,error_message); */ /* class_alloc(del,sizeof(double)*neqp,error_message); */ /* class_alloc(yinterp,sizeof(double)*neqp,error_message); */ /* class_alloc(ypinterp,sizeof(double)*neqp,error_message); */ /* class_alloc(yppinterp,sizeof(double)*neqp,error_message); */ /* class_alloc(tempvec1,sizeof(double)*neqp,error_message); */ /* class_alloc(tempvec2,sizeof(double)*neqp,error_message); */ /* class_alloc(interpidx,sizeof(int)*neqp,error_message); */ /* Allocate vector of pointers to rows of dif:*/ /* class_alloc(dif,sizeof(double*)*neqp,error_message); */ /* class_calloc(dif[1],(7*neq+1),sizeof(double),error_message); */ /* dif[0] = NULL; */ /* for(j=2;j<=neq;j++) dif[j] = dif[j-1]+7; */ /* Set row pointers... */ /*Set pointers:*/ ynew = y_inout-1; /* This way y_inout is always up to date. */ /*Initialize the jacobian:*/ class_call(initialize_jacobian(&jac,neq,error_message),error_message,error_message); /* Initialize workspace for numjac: */ class_call(initialize_numjac_workspace(&nj_ws,neq,error_message),error_message,error_message); /* Initialize some method parameters:*/ for(ii=0;ii<5;ii++){ invGa[ii] = 1.0/(G[ii]*(1.0 - alpha[ii])); erconst[ii] = alpha[ii]*G[ii] + 1.0/(2.0+ii); } /* Set the relevant indices which needs to be found by interpolation. */ /* But if we want to print variables for testing purposes, just interpolate everything.. */ for(ii=1;ii<=neq;ii++){ y[ii] = y_inout[ii-1]; if (print_variables==NULL){ interpidx[ii]=used_in_output[ii-1]; } else{ interpidx[ii]=1; } } t0 = x_ini; tfinal = x_final; /* Some CLASS-specific stuff:*/ next=0; while (t_vec[next] < t0) next++; if (verbose > 1){ numidx=0; for(ii=1;ii<=neq;ii++){ if (interpidx[ii]==_TRUE_) numidx++; } printf("%d/%d ",numidx,neq); } htspan = fabs(tfinal-t0); for(ii=0;ii<6;ii++) stepstat[ii] = 0; class_call((*derivs)(t0,y+1,f0+1,parameters_and_workspace_for_derivs,error_message),error_message,error_message); stepstat[2] +=1; if ((tfinal-t0)<0.0){ tdir = -1; } else{ tdir = 1; } hmax = (tfinal-t0)/10.0; t = t0; nfenj=0; class_call(numjac((*derivs),t,y,f0,&jac,&nj_ws,abstol,neq, &nfenj,parameters_and_workspace_for_derivs,error_message), error_message,error_message); stepstat[3] += 1; stepstat[2] += nfenj; Jcurrent = _TRUE_; /* True */ hmin = 1.e-20;//16.0*eps*fabs(t); /*Calculate initial step */ rh = 0.0; for(jj=1;jj<=neq;jj++){ wt[jj] = MAX(fabs(y[jj]),threshold); /*printf("wt: %4.8f \n",wt[jj]);*/ rh = MAX(rh,1.25/sqrt(rtol)*fabs(f0[jj]/wt[jj])); } absh = MIN(hmax, htspan); if (absh * rh > 1.0) absh = 1.0 / rh; absh = MAX(absh, hmin); h = tdir * absh; tdel = (t + tdir*MIN(sqrt(eps)*MAX(fabs(t),fabs(t+h)),absh)) - t; class_call((*derivs)(t+tdel,y+1,tempvec1+1,parameters_and_workspace_for_derivs,error_message), error_message,error_message); stepstat[2] += 1; /*I assume that a full jacobi matrix is always calculated in the beginning...*/ for(ii=1;ii<=neq;ii++){ ddfddt[ii]=0.0; for(jj=1;jj<=neq;jj++){ ddfddt[ii]+=(jac.dfdy[ii][jj])*f0[jj]; } } rh = 0.0; for(ii=1;ii<=neq;ii++){ ddfddt[ii] += (tempvec1[ii] - f0[ii]) / tdel; rh = MAX(rh,1.25*sqrt(0.5*fabs(ddfddt[ii]/wt[ii])/rtol)); } absh = MIN(hmax, htspan); if (absh * rh > 1.0) absh = 1.0 / rh; absh = MAX(absh, hmin); h = tdir * absh; /* Done calculating initial step Get ready to do the loop:*/ k = 1; /*start at order 1 with BDF1 */ klast = k; abshlast = absh; for(ii=1;ii<=neq;ii++) dif[ii][1] = h*f0[ii]; hinvGak = h*invGa[k-1]; nconhk = 0; /*steps taken with current h and k*/ class_call(new_linearisation(&jac,hinvGak,neq,error_message), error_message,error_message); stepstat[4] += 1; havrate = _FALSE_; /*false*/ /* Doing main loop: */ done = _FALSE_; at_hmin = _FALSE_; while (done==_FALSE_){ hmin = minimum_variation; maxtmp = MAX(hmin,absh); absh = MIN(hmax, maxtmp); if (fabs(absh-hmin)<100*eps){ /* If the stepsize has not changed */ if (at_hmin==_TRUE_){ absh = abshlast; /*required by stepsize recovery */ } at_hmin = _TRUE_; } else{ at_hmin = _FALSE_; } h = tdir * absh; /* Stretch the step if within 10% of tfinal-t. */ if (1.1*absh >= fabs(tfinal - t)){ h = tfinal - t; absh = fabs(h); done = _TRUE_; } if (((fabs(absh-abshlast)/absh)>1e-6)||(k!=klast)){ adjust_stepsize(dif,(absh/abshlast),neq,k); hinvGak = h * invGa[k-1]; nconhk = 0; class_call(new_linearisation(&jac,hinvGak,neq,error_message), error_message,error_message); stepstat[4] += 1; havrate = _FALSE_; } /* Loop for advancing one step */ nofailed = _TRUE_; for( ; ; ){ gotynew = _FALSE_; /* is ynew evaluated yet?*/ while(gotynew==_FALSE_){ /*Compute the constant terms in the equation for ynew. Next FOR lop is just: psi = matmul(dif(:,1:k),(G(1:k) * invGa(k)))*/ for(ii=1;ii<=neq;ii++){ psi[ii] = 0.0; for(jj=1;jj<=k;jj++){ psi[ii] += dif[ii][jj]*G[jj-1]*invGa[k-1]; } } /* Predict a solution at t+h. */ tnew = t + h; if (done==_TRUE_){ tnew = tfinal; /*Hit end point exactly. */ } h = tnew - t; /* Purify h. */ for(ii=1;ii<=neq;ii++){ pred[ii] = y[ii]; for(jj=1;jj<=k;jj++){ pred[ii] +=dif[ii][jj]; } } eqvec(pred,ynew,neq); /*The difference, difkp1, between pred and the final accepted ynew is equal to the backward difference of ynew of order k+1. Initialize to zero for the iteration to compute ynew. */ minnrm = 0.0; for(j=1;j<=neq;j++){ difkp1[j] = 0.0; maxtmp = MAX(fabs(ynew[j]),fabs(y[j])); invwt[j] = 1.0 / MAX(maxtmp,threshold); maxtmp = 100*eps*fabs(ynew[j]*invwt[j]); minnrm = MAX(minnrm,maxtmp); } /* Iterate with simplified Newton method. */ tooslow = _FALSE_; for(iter=1;iter<=maxit;iter++){ for (ii=1;ii<=neq;ii++){ tempvec1[ii]=(psi[ii]+difkp1[ii]); } class_call((*derivs)(tnew,ynew+1,f0+1,parameters_and_workspace_for_derivs,error_message), error_message,error_message); stepstat[2] += 1; for(j=1;j<=neq;j++){ rhs[j] = hinvGak*f0[j]-tempvec1[j]; } /*Solve the linear system A*x=del by using the LU decomposition stored in jac.*/ if (jac.use_sparse){ sp_lusolve(jac.Numerical, rhs+1, del+1); } else{ eqvec(rhs,del,neq); lubksb(jac.LU,neq,jac.luidx,del); } stepstat[5]+=1; newnrm = 0.0; for(j=1;j<=neq;j++){ maxtmp = fabs(del[j]*invwt[j]); newnrm = MAX(newnrm,maxtmp); } for(j=1;j<=neq;j++){ difkp1[j] += del[j]; ynew[j] = pred[j] + difkp1[j]; } if (newnrm <= minnrm){ gotynew = _TRUE_; break; /* Break Newton loop */ } else if(iter == 1){ if (havrate==_TRUE_){ errit = newnrm * rate / (1.0 - rate); if (errit <= 0.05*rtol){ gotynew = _TRUE_; break; /* Break Newton Loop*/ } } else { rate = 0.0; } } else if(newnrm > 0.9*oldnrm){ tooslow = _TRUE_; break; /*Break Newton lop */ } else{ rate = MAX(0.9*rate, newnrm / oldnrm); havrate = _TRUE_; errit = newnrm * rate / (1.0 - rate); if (errit <= 0.5*rtol){ gotynew = _TRUE_; break; /* exit newton */ } else if (iter == maxit){ tooslow = _TRUE_; break; /*exit newton */ } else if (0.5*rtol < errit*pow(rate,(maxit-iter))){ tooslow = _TRUE_; break; /*exit Newton */ } } oldnrm = newnrm; } if (tooslow==_TRUE_){ stepstat[1] += 1; /* ! Speed up the iteration by forming new linearization or reducing h. */ if (Jcurrent==_FALSE_){ class_call((*derivs)(t,y+1,f0+1,parameters_and_workspace_for_derivs,error_message), error_message,error_message); nfenj=0; class_call(numjac((*derivs),t,y,f0,&jac,&nj_ws,abstol,neq, &nfenj,parameters_and_workspace_for_derivs,error_message), error_message,error_message); stepstat[3] += 1; stepstat[2] += (nfenj + 1); Jcurrent = _TRUE_; } else if (absh <= hmin){ class_test(absh <= hmin, error_message, "Step size too small: step:%g, minimum:%g, in interval: [%g:%g]\n", absh,hmin,t0,tfinal); } else{ abshlast = absh; absh = MAX(0.3 * absh, hmin); h = tdir * absh; done = _FALSE_; adjust_stepsize(dif,(absh/abshlast),neq,k); hinvGak = h * invGa[k-1]; nconhk = 0; } /* A new linearisation is needed in both cases */ class_call(new_linearisation(&jac,hinvGak,neq,error_message), error_message,error_message); stepstat[4] += 1; havrate = _FALSE_; } } /*end of while loop for getting ynew difkp1 is now the backward difference of ynew of order k+1. */ err = 0.0; for(jj=1;jj<=neq;jj++){ err = MAX(err,fabs(difkp1[jj]*invwt[jj])); } err = err * erconst[k-1]; if (err>rtol){ /*Step failed */ stepstat[1]+= 1; if (absh <= hmin){ //BEN FLAG: I REMOVED THIS FOR NO GOOD REASON!/ class_test(absh <= hmin, error_message, "Step size too small: step:%g, minimum:%g, in interval: [%g:%g]\n", absh,hmin,t0,tfinal); } abshlast = absh; if (nofailed==_TRUE_){ nofailed = _FALSE_; hopt = absh * MAX(0.1, 0.833*pow((rtol/err),(1.0/(k+1)))); if (k > 1){ errkm1 = 0.0; for(jj=1;jj<=neq;jj++){ errkm1 = MAX(errkm1,fabs((dif[jj][k]+difkp1[jj])*invwt[jj])); } errkm1 = errkm1*erconst[k-2]; hkm1 = absh * MAX(0.1, 0.769*pow((rtol/errkm1),(1.0/k))); if (hkm1 > hopt){ hopt = MIN(absh,hkm1); /* don't allow step size increase */ k = k - 1; } } absh = MAX(hmin, hopt); } else{ absh = MAX(hmin, 0.5 * absh); } h = tdir * absh; if (absh < abshlast){ done = _FALSE_; } adjust_stepsize(dif,(absh/abshlast),neq,k); hinvGak = h * invGa[k-1]; nconhk = 0; class_call(new_linearisation(&jac,hinvGak,neq,error_message), error_message,error_message); stepstat[4] += 1; havrate = _FALSE_; } else { break; /* Succesfull step */ } } /* End of conditionless FOR loop */ stepstat[0] += 1; /* Update dif: */ for(jj=1;jj<=neq;jj++){ dif[jj][k+2] = difkp1[jj] - dif[jj][k+1]; dif[jj][k+1] = difkp1[jj]; } for(j=k;j>=1;j--){ for(ii=1;ii<=neq;ii++){ dif[ii][j] += dif[ii][j+1]; } } /** Output **/ while ((next<tres)&&(tdir * (tnew - t_vec[next]) >= 0.0)){ /* Do we need to write output? */ if (tnew==t_vec[next]){ class_call((*output)(t_vec[next],ynew+1,f0+1,next,parameters_and_workspace_for_derivs,error_message), error_message,error_message); // MODIFICATION BY LUC // All print_variables have been moved to the end of time step /* if (print_variables != NULL){ class_call((*print_variables)(t_vec[next],ynew+1,f0+1, parameters_and_workspace_for_derivs,error_message), error_message,error_message); } */ } else { /*Interpolate if we have overshot sample values*/ interp_from_dif(t_vec[next],tnew,ynew,h,dif,k,yinterp,ypinterp,yppinterp,interpidx,neq,2); class_call((*output)(t_vec[next],yinterp+1,ypinterp+1,next,parameters_and_workspace_for_derivs, error_message),error_message,error_message); } next++; } /** End of output **/ if (done==_TRUE_) { break; } klast = k; abshlast = absh; nconhk = MIN(nconhk+1,maxk+2); if (nconhk >= k + 2){ temp = 1.2*pow((err/rtol),(1.0/(k+1.0))); if (temp > 0.1){ hopt = absh / temp; } else { hopt = 10*absh; } kopt = k; if (k > 1){ errkm1 = 0.0; for(jj=1;jj<=neq;jj++){ errkm1 = MAX(errkm1,fabs(dif[jj][k]*invwt[jj])); } errkm1 = errkm1*erconst[k-2]; temp = 1.3*pow((errkm1/rtol),(1.0/k)); if (temp > 0.1){ hkm1 = absh / temp; } else { hkm1 = 10*absh; } if (hkm1 > hopt){ hopt = hkm1; kopt = k - 1; } } if (k < maxk){ errkp1 = 0.0; for(jj=1;jj<=neq;jj++){ errkp1 = MAX(errkp1,fabs(dif[jj][k+2]*invwt[jj])); } errkp1 = errkp1*erconst[k]; temp = 1.4*pow((errkp1/rtol),(1.0/(k+2.0))); if (temp > 0.1){ hkp1 = absh / temp; } else { hkp1 = 10*absh; } if (hkp1 > hopt){ hopt = hkp1; kopt = k + 1; } } if (hopt > absh){ absh = hopt; if (k!=kopt){ k = kopt; } } } /* Advance the integration one step. */ t = tnew; eqvec(ynew,y,neq); Jcurrent = _FALSE_; // MODIFICATION BY LUC if (print_variables!=NULL){ class_call((*derivs)(tnew, ynew+1, f0+1, parameters_and_workspace_for_derivs,error_message), error_message, error_message); class_call((*print_variables)(tnew,ynew+1,f0+1, parameters_and_workspace_for_derivs,error_message), error_message,error_message); } // end of modification } /* a last call is compulsory to ensure that all quantitites in y,dy,parameters_and_workspace_for_derivs are updated to the last point in the covered range */ class_call( (*derivs)(tnew, ynew+1, f0+1, parameters_and_workspace_for_derivs,error_message), error_message, error_message); if (verbose > 0){ printf("\n End of evolver. Next=%d, t=%e and tnew=%e.",next,t,tnew); printf("\n Statistics: [%d %d %d %d %d %d] \n",stepstat[0],stepstat[1], stepstat[2],stepstat[3],stepstat[4],stepstat[5]); } /** Deallocate memory */ free(buffer); /* free(f0); */ /* free(wt); */ /* free(ddfddt); */ /* free(pred); */ /* free(y); */ /* free(invwt); */ /* free(rhs); */ /* free(psi); */ /* free(difkp1); */ /* free(del); */ /* free(yinterp); */ /* free(ypinterp); */ /* free(yppinterp); */ /* free(tempvec1); */ /* free(tempvec2); */ /* free(interpidx); */ /* free(dif[1]); */ /* free(dif); */ uninitialize_jacobian(&jac); uninitialize_numjac_workspace(&nj_ws); return _SUCCESS_; } /*End of program*/