/*! \fn void MultiNSH(int n, Real *tstop, Real *mratio, Real etavk, * Real *uxNSH, Real *uyNSH, Real *wxNSH, Real *wyNSH) * \brief Multi-species NSH equilibrium * * Input: # of particle types (n), dust stopping time and mass ratio array, and * drift speed etavk. * Output: gas NSH equlibrium velocity (u), and dust NSH equilibrium velocity * array (w). */ void MultiNSH(int n, Real *tstop, Real *mratio, Real etavk, Real *uxNSH, Real *uyNSH, Real *wxNSH, Real *wyNSH) { int i,j; Real *Lambda1,**Lam1GamP1, **A, **B, **Tmp; Lambda1 = (Real*)calloc_1d_array(n, sizeof(Real)); /* Lambda^{-1} */ Lam1GamP1=(Real**)calloc_2d_array(n, n, sizeof(Real)); /* Lambda1*(1+Gamma) */ A = (Real**)calloc_2d_array(n, n, sizeof(Real)); B = (Real**)calloc_2d_array(n, n, sizeof(Real)); Tmp = (Real**)calloc_2d_array(n, n, sizeof(Real)); /* definitions */ for (i=0; i<n; i++){ for (j=0; j<n; j++) Lam1GamP1[i][j] = mratio[j]; Lam1GamP1[i][i] += 1.0; Lambda1[i] = 1.0/(tstop[i]+1.0e-16); for (j=0; j<n; j++) Lam1GamP1[i][j] *= Lambda1[i]; } /* Calculate A and B */ MatrixMult(Lam1GamP1, Lam1GamP1, n,n,n, Tmp); for (i=0; i<n; i++) Tmp[i][i] += 1.0; InverseMatrix(Tmp, n, B); for (i=0; i<n; i++) for (j=0; j<n; j++) B[i][j] *= Lambda1[j]; MatrixMult(Lam1GamP1, B, n,n,n, A); /* Obtain NSH velocities */ *uxNSH = 0.0; *uyNSH = 0.0; for (i=0; i<n; i++){ wxNSH[i] = 0.0; wyNSH[i] = 0.0; for (j=0; j<n; j++){ wxNSH[i] -= B[i][j]; wyNSH[i] -= A[i][j]; } wxNSH[i] *= 2.0*etavk; wyNSH[i] *= etavk; *uxNSH -= mratio[i]*wxNSH[i]; *uyNSH -= mratio[i]*wyNSH[i]; wyNSH[i] += etavk; } free(Lambda1); free_2d_array(A); free_2d_array(B); free_2d_array(Lam1GamP1); free_2d_array(Tmp); return; }
/*! \fn void integrate_init_2d(MeshS *pM) * \brief Allocate temporary integration arrays */ void integrate_init_2d(MeshS *pM) { int nmax,size1=0,size2=0,size3=0,nl,nd; /* Cycle over all Grids on this processor to find maximum Nx1, Nx2, Nx3 */ for (nl=0; nl<(pM->NLevels); nl++){ for (nd=0; nd<(pM->DomainsPerLevel[nl]); nd++){ if (pM->Domain[nl][nd].Grid != NULL) { if (pM->Domain[nl][nd].Grid->Nx[0] > size1){ size1 = pM->Domain[nl][nd].Grid->Nx[0]; } if (pM->Domain[nl][nd].Grid->Nx[1] > size2){ size2 = pM->Domain[nl][nd].Grid->Nx[1]; } if (pM->Domain[nl][nd].Grid->Nx[2] > size3){ size3 = pM->Domain[nl][nd].Grid->Nx[2]; } } } } size1 = size1 + 2*nghost; size2 = size2 + 2*nghost; size3 = size3 + 2*nghost; nmax = MAX((MAX(size1,size2)),size3); /*refer to material integrate_2d_ctu.c*/ if ((Bxc = (Real*)malloc(nmax*sizeof(Real))) == NULL) goto on_error; if ((Bxi = (Real*)malloc(nmax*sizeof(Real))) == NULL) goto on_error; if ((U1d= (Cons1DS*)malloc(nmax*sizeof(Cons1DS))) == NULL) goto on_error; if ((W = (Prim1DS*)malloc(nmax*sizeof(Prim1DS))) == NULL) goto on_error; if ((x1Flux =(Cons1DS**)calloc_2d_array(size2,size1,sizeof(Cons1DS)))==NULL) goto on_error; if ((x2Flux =(Cons1DS**)calloc_2d_array(size2,size1,sizeof(Cons1DS)))==NULL) goto on_error; return; on_error: integrate_destruct(); ath_error("[integrate_init]: malloc returned a NULL pointer\n"); }
/*! \fn static void flux_output_func(MeshS *pM, OutputS *pOut) * \brief New output format which outputs y-integrated angular momentum fluxes * Currently can only be used with 1 proc and 1 domain. */ static void flux_output_func(MeshS *pM, OutputS *pOut) { GridS *pG=pM->Domain[0][0].Grid; int nx1,nx2,nx3, ind, i,j,k, ind_i,ind_k; Real lx1, lx2, lx3; PrimS W[7],Ws[7],We[7]; Real x1[7],x2[7],x3[7],xs1[7],xs2[7],xs3[7],xe1[7],xe2[7],xe3[7]; Real dmin, dmax; Real **Fluxx=NULL; Real **FluxH=NULL; Real **FluxNu=NULL; Real **Th=NULL; Real **outCoordsx1=NULL; Real **outCoordsx3=NULL; Real **vx=NULL; Real **vy=NULL; Real **sigvx=NULL; Real **osigx=NULL; Real **davg=NULL; FILE *pfile; char *fname; nx1 = pG->Nx[0]; nx3 = pG->Nx[2]; nx2 = pG->Nx[1]; lx1 = pG->MaxX[0] - pG->MinX[0]; lx2 = pG->MaxX[1] - pG->MinX[1]; lx3 = pG->MaxX[2] - pG->MinX[2]; printf("%d, %d, %d, %g, %g, %g\n",nx1,nx2,nx3,lx1,lx2,lx3); #ifdef MPI_PARALLEL printf("%d, %d, %d, %g, %g, %g \n", myID_Comm_world,nx1,nx3,lx1,lx2,lx3); printf("IND %d: (%d,%d), (%d,%d)\n",myID_Comm_world,pG->is,pG->js,pG->ie,pG->je); #endif Fluxx=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (Fluxx == NULL) return; FluxH=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (FluxH == NULL) return; FluxNu=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (FluxNu == NULL) return; Th=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (Th == NULL) return; outCoordsx1=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (outCoordsx1 == NULL) return; outCoordsx3=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (outCoordsx3 == NULL) return; vx=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (vx == NULL) return; vy=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (vy == NULL) return; sigvx=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (sigvx == NULL) return; osigx=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (osigx == NULL) return; davg=(Real **)calloc_2d_array(nx3,nx1,sizeof(Real)); if (davg == NULL) return; /* Open file and write header */ if((fname = ath_fname(NULL,pM->outfilename,NULL,NULL,num_digit, pOut->num,pOut->id,"tab")) == NULL){ ath_error("[dump_tab]: Error constructing filename\n"); } if((pfile = fopen(fname,"w")) == NULL){ ath_error("[dump_tab]: Unable to open ppm file %s\n",fname); } free(fname); #ifdef MPI_PARALLEL if(myID_Comm_world == 0) #endif fprintf(pfile,"#t=%12.8e x,FH,Fx,Fnu,Th,vx,vy,davg,sigvx,omsigx \n", pOut->t); /* Compute y-integrated fluxes explicitly. * For the derivatives use a central difference method. * For the integration use a composite trapezoid method. * Both of these make use of the ghost cells for the boundary values. * 1 FluxH = < d * vx1 *vx2 > * 2 Fluxx = < 0.5 * omega * d * x1 * vx1 > * 3 FluxNu = - < nu_iso * d/dx vx2 > * 4 Th = - < d * d/dy phi > * 5 vx = < vx1 > * 6 vy = < vx2 > * 7 davg = < d > * 8 sigvx = < d * vx1 > * 9 osigx = < 0.5 * sig * omega * x1 > * * * x 4 x * 1 0 2 * x 3 x * * The variables are stored up as arrays of primitives W[7]. * W = ( W[k][j][i], W[k][j][i-1], W[k][j][i+1], W[k][j-1][i], * W[k][j+1][i], W[k-1][j][i], W[k+1][j][i] ) * This is needed for the integration and differentiation. * * The background shear is taken out of vy when doing computations. */ for(k=pG->ks; k <=pG->ke; k++) { for(i=pG->is; i<= pG->ie; i++) { ind_i=i-pG->is; ind_k=k-pG->ks; for(ind=0; ind < 7; ind++) { if (ind==0) { cc_pos(pG,i,pG->js,k,&xs1[ind],&xs2[ind],&xs3[ind]); cc_pos(pG,i,pG->je,k,&xe1[ind],&xe2[ind],&xe3[ind]); Ws[ind] = Cons_to_Prim(&(pG->U[k][pG->js][i])); We[ind] = Cons_to_Prim(&(pG->U[k][pG->je][i])); } if (ind > 0 && ind < 3) { cc_pos(pG,i+2*ind-3,pG->js,k,&xs1[ind],&xs2[ind],&xs3[ind]); cc_pos(pG,i+2*ind-3,pG->je,k,&xe1[ind],&xe2[ind],&xe3[ind]); Ws[ind] = Cons_to_Prim(&(pG->U[k][pG->js][i+2*ind-3])); We[ind] = Cons_to_Prim(&(pG->U[k][pG->je][i+2*ind-3])); } if (ind > 2 && ind < 5) { cc_pos(pG,i,pG->js+2*ind-7,k,&xs1[ind],&xs2[ind],&xs3[ind]); cc_pos(pG,i,pG->je+2*ind-7,k,&xe1[ind],&xe2[ind],&xe3[ind]); Ws[ind] = Cons_to_Prim(&(pG->U[k][pG->js+2*ind-7][i])); We[ind] = Cons_to_Prim(&(pG->U[k][pG->je+2*ind-7][i])); } if (ind > 4 && ind < 7) { if (pG->MinX[2] == pG->MaxX[2]) { cc_pos(pG,i,pG->js,k,&xs1[ind],&xs2[ind],&xs3[ind]); cc_pos(pG,i,pG->je,k,&xe1[ind],&xe2[ind],&xe3[ind]); Ws[ind] = Cons_to_Prim(&(pG->U[k][pG->js][i])); We[ind] = Cons_to_Prim(&(pG->U[k][pG->je][i])); } else { cc_pos(pG,i,pG->js,k+2*ind-11,&xs1[ind],&xs2[ind],&xs3[ind]); cc_pos(pG,i,pG->je,k+2*ind-11,&xe1[ind],&xe2[ind],&xe3[ind]); Ws[ind] = Cons_to_Prim(&(pG->U[k+2*ind-11][pG->js][i])); We[ind] = Cons_to_Prim(&(pG->U[k+2*ind-11][pG->je][i])); } } Ws[ind].V2 = Ws[ind].V2 + qshear*Omega_0*xs1[ind]; We[ind].V2 = We[ind].V2 + qshear*Omega_0*xe1[ind]; /* If using d' instead of d then put in * W[ind] = W[ind]->d - d0; */ } /* Set initial values for the integration */ FluxH[ind_k][ind_i] = 0.5*( Ws[0].d * Ws[0].V1 * Ws[0].V2 + We[0].d * We[0].V1 * We[0].V2 ); Fluxx[ind_k][ind_i] = 0.5*( Ws[0].d * Ws[0].V1 * xs1[0] + We[0].d * We[0].V1 * xe1[0] ); #ifdef VISCOSITY FluxNu[ind_k][ind_i] = 0.5*( Ws[2].V2 - Ws[1].V2 + We[2].V2 - We[1].V2 ); #else FluxNu[ind_k][ind_i] = 0; #endif Th[ind_k][ind_i] = 0.5*( Ws[0].d * dx2PlanetPot(xs1[0],xs2[0],xs3[0]) + We[0].d * dx2PlanetPot(xe1[0],xe2[0],xe3[0]) ); vx[ind_k][ind_i] = 0.5*( Ws[0].V1 + We[0].V1 ); vy[ind_k][ind_i] = 0.5*( Ws[0].V2 + We[0].V2 ); sigvx[ind_k][ind_i] = 0.5*( Ws[0].d * Ws[0].V1 + We[0].d * We[0].V1 ); osigx[ind_k][ind_i] = 0.5*( Ws[0].d * xs1[0] + We[0].d * xe1[0] ); davg[ind_k][ind_i] = 0.5*(Ws[0].d + We[0].d); for(j=(pG->js)+1; j < pG->je; j++) { for(ind=0; ind < 7; ind++) { if (ind==0) { cc_pos(pG,i,j,k,&x1[ind],&x2[ind],&x3[ind]); W[ind] = Cons_to_Prim(&(pG->U[k][j][i])); } if (ind > 0 && ind < 3) { cc_pos(pG,i+2*ind-3,j,k,&x1[ind],&x2[ind],&x3[ind]); W[ind] = Cons_to_Prim(&(pG->U[k][j][i+2*ind-3])); } if (ind > 2 && ind < 5) { cc_pos(pG,i,j+2*ind-7,k,&x1[ind],&x2[ind],&x3[ind]); W[ind] = Cons_to_Prim(&(pG->U[k][j+2*ind-7][i])); } if (ind > 4 && ind < 7) { if (pG->MinX[2] == pG->MaxX[2]) { cc_pos(pG,i,j,k,&x1[ind],&x2[ind],&x3[ind]); W[ind] = Cons_to_Prim(&(pG->U[k][j][i])); } else { cc_pos(pG,i,j,k+2*ind-11,&x1[ind],&x2[ind],&x3[ind]); W[ind] = Cons_to_Prim(&(pG->U[k+2*ind-11][j][i])); } } W[ind].V2 = W[ind].V2 + qshear*Omega_0*x1[ind]; } FluxH[ind_k][ind_i] += W[0].d * W[0].V1 * W[0].V2; Fluxx[ind_k][ind_i] += W[0].d * W[0].V1 * x1[0]; #ifdef VISCOSITY FluxNu[ind_k][ind_i] += W[0].d * (W[2].V2 - W[1].V1); #endif Th[ind_k][ind_i] += W[0].d * dx2PlanetPot(x1[0],x2[0],x3[0]); vx[ind_k][ind_i] += W[0].V1; vy[ind_k][ind_i] += W[0].V2; sigvx[ind_k][ind_i] += W[0].d * W[0].V1; osigx[ind_k][ind_i] += W[0].d * x1[0]; davg[ind_k][ind_i] += W[0].d; } FluxH[ind_k][ind_i] *= (pG->dx2)/lx2; Fluxx[ind_k][ind_i] *= .5*Omega_0*(pG->dx2)/lx2; #ifdef VISCOSITY FluxNu[ind_k][ind_i] *= -nu_iso*(pG->dx2)/(2*lx2*(pG->dx1)); #endif Th[ind_k][ind_i] *= -1.0*(pG->dx2)/lx2; vx[ind_k][ind_i] *= (pG->dx2)/lx2; vy[ind_k][ind_i] *= (pG->dx2)/lx2; sigvx[ind_k][ind_i] *= (pG->dx2)/lx2; osigx[ind_k][ind_i] *= (0.5*Omega_0*(pG->dx2))/lx2; davg[ind_k][ind_i] *= (pG->dx2)/lx2; outCoordsx1[ind_k][ind_i]=x1[0]; outCoordsx3[ind_k][ind_i]=x3[0]; } } /* Quantities are ready to be written to output file * Format (Not outputting x3 at the moment: * x1 (x3) FluxH Fluxx Fluxnu Th vx vy davg sigvx osigx */ for(k=pG->ks; k<=pG->ke; k++) { for(i=pG->is; i<=pG->ie; i++) { ind_k=k-pG->ks; ind_i=i-pG->is; if (lx3==0) { fprintf(pfile,"%12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e\n", outCoordsx1[ind_k][ind_i],FluxH[ind_k][ind_i], Fluxx[ind_k][ind_i],FluxNu[ind_k][ind_i],Th[ind_k][ind_i], vx[ind_k][ind_i],vy[ind_k][ind_i],davg[ind_k][ind_i],sigvx[ind_k][ind_i], osigx[ind_k][ind_i]); } else { fprintf(pfile,"%12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e %12.8e\n", outCoordsx1[ind_k][ind_i],outCoordsx3[ind_k][ind_i],FluxH[ind_k][ind_i], Fluxx[ind_k][ind_i],FluxNu[ind_k][ind_i],Th[ind_k][ind_i], vx[ind_k][ind_i],vy[ind_k][ind_i],davg[ind_k][ind_i],sigvx[ind_k][ind_i], osigx[ind_k][ind_i]); } } } fclose(pfile); free_2d_array(Fluxx); free_2d_array(FluxH); free_2d_array(FluxNu); free_2d_array(Th); free_2d_array(outCoordsx1); free_2d_array(outCoordsx3); free_2d_array(vx); free_2d_array(vy); free_2d_array(sigvx); free_2d_array(osigx); free_2d_array(davg); return; }
void problem(DomainS *pDomain) { GridS *pGrid = pDomain->Grid; ConsS **Soln; int i, is = pGrid->is, ie = pGrid->ie; int j, js = pGrid->js, je = pGrid->je; int k, ks = pGrid->ks, ke = pGrid->ke; int nx1, nx2; int dir; Real angle; /* Angle the wave direction makes with the x1-direction */ Real x1size,x2size,x1,x2,x3,cs,sn; Real v_par, v_perp, den, pres; Real lambda; /* Wavelength */ #ifdef RESISTIVITY Real v_A, kva, omega_h, omega_l, omega_r; #endif nx1 = (ie - is + 1) + 2*nghost; nx2 = (je - js + 1) + 2*nghost; if (pGrid->Nx[1] == 1) { ath_error("[problem] Grid must be 2D"); } if ((Soln = (ConsS**)calloc_2d_array(nx2,nx1,sizeof(ConsS))) == NULL) ath_error("[problem]: Error allocating memory for Soln\n"); if (pDomain->Level == 0){ if ((RootSoln =(ConsS**)calloc_2d_array(nx2,nx1,sizeof(ConsS)))==NULL) ath_error("[problem]: Error allocating memory for RootSoln\n"); } /* An angle = 0.0 is a wave aligned with the x1-direction. */ /* An angle = 90.0 is a wave aligned with the x2-direction. */ angle = par_getd("problem","angle"); x1size = pDomain->RootMaxX[0] - pDomain->RootMinX[0]; x2size = pDomain->RootMaxX[1] - pDomain->RootMinX[1]; /* Compute the sin and cos of the angle and the wavelength. */ /* Put one wavelength in the grid */ if (angle == 0.0) { sin_a = 0.0; cos_a = 1.0; lambda = x1size; } else if (angle == 90.0) { sin_a = 1.0; cos_a = 0.0; lambda = x2size; } else { /* We put 1 wavelength in each direction. Hence the wavelength * lambda = (pDomain->RootMaxX[0] - pDomain->RootMinX[0])*cos_a * AND lambda = (pDomain->RootMaxX[1] - pDomain->RootMinX[1])*sin_a; * are both satisfied. */ if(x1size == x2size){ cos_a = sin_a = sqrt(0.5); } else{ angle = atan((double)(x1size/x2size)); sin_a = sin(angle); cos_a = cos(angle); } /* Use the larger angle to determine the wavelength */ if (cos_a >= sin_a) { lambda = x1size*cos_a; } else { lambda = x2size*sin_a; } } /* Initialize k_parallel */ k_par = 2.0*PI/lambda; b_par = par_getd("problem","b_par"); den = 1.0; ath_pout(0,"va_parallel = %g\nlambda = %g\n",b_par/sqrt(den),lambda); b_perp = par_getd("problem","b_perp"); v_perp = b_perp/sqrt((double)den); dir = par_geti_def("problem","dir",1); /* right(1)/left(2) polarization */ if (dir == 1) /* right polarization */ fac = 1.0; else /* left polarization */ fac = -1.0; #ifdef RESISTIVITY Q_Hall = par_getd("problem","Q_H"); d_ind = 0.0; v_A = b_par/sqrt((double)den); if (Q_Hall > 0.0) { kva = k_par*v_A; omega_h = 1.0/Q_Hall; omega_r = 0.5*SQR(kva)/omega_h*(sqrt(1.0+SQR(2.0*omega_h/kva)) + 1.0); omega_l = 0.5*SQR(kva)/omega_h*(sqrt(1.0+SQR(2.0*omega_h/kva)) - 1.0); if (dir == 1) /* right polarization (whistler wave) */ v_perp = v_perp * kva / omega_r; else /* left polarization */ v_perp = v_perp * kva / omega_l; } #endif /* The gas pressure and parallel velocity are free parameters. */ pres = par_getd("problem","pres"); v_par = par_getd("problem","v_par"); /* Use the vector potential to initialize the interface magnetic fields * The iterface fields are located at the left grid cell face normal */ for (k=ks; k<=ke; k++) { for (j=js; j<=je+1; j++) { for (i=is; i<=ie+1; i++) { cc_pos(pGrid,i,j,k,&x1,&x2,&x3); cs = cos(k_par*(x1*cos_a + x2*sin_a)); x1 -= 0.5*pGrid->dx1; x2 -= 0.5*pGrid->dx2; pGrid->B1i[k][j][i] = -(A3(x1,(x2+pGrid->dx2)) - A3(x1,x2))/pGrid->dx2; pGrid->B2i[k][j][i] = (A3((x1+pGrid->dx1),x2) - A3(x1,x2))/pGrid->dx1; pGrid->B3i[k][j][i] = b_perp*cs; } } } if (pGrid->Nx[2] > 1) { for (j=js; j<=je+1; j++) { for (i=is; i<=ie+1; i++) { cc_pos(pGrid,i,j,k,&x1,&x2,&x3); cs = cos(k_par*(x1*cos_a + x2*sin_a)); pGrid->B3i[ke+1][j][i] = b_perp*cs; } } } /* Now initialize the cell centered quantities */ for (k=ks; k<=ke; k++) { for (j=js; j<=je; j++) { for (i=is; i<=ie; i++) { cc_pos(pGrid,i,j,k,&x1,&x2,&x3); sn = sin(k_par*(x1*cos_a + x2*sin_a)) * fac; cs = cos(k_par*(x1*cos_a + x2*sin_a)); Soln[j][i].d = den; Soln[j][i].M1 = den*(v_par*cos_a + v_perp*sn*sin_a); Soln[j][i].M2 = den*(v_par*sin_a - v_perp*sn*cos_a); Soln[j][i].M3 = -den*v_perp*cs; pGrid->U[k][j][i].d = Soln[j][i].d; pGrid->U[k][j][i].M1 = Soln[j][i].M1; pGrid->U[k][j][i].M2 = Soln[j][i].M2; pGrid->U[k][j][i].M3 = Soln[j][i].M3; Soln[j][i].B1c = 0.5*(pGrid->B1i[k][j][i] + pGrid->B1i[k][j][i+1]); Soln[j][i].B2c = 0.5*(pGrid->B2i[k][j][i] + pGrid->B2i[k][j+1][i]); Soln[j][i].B3c = b_perp*cs; pGrid->U[k][j][i].B1c = Soln[j][i].B1c; pGrid->U[k][j][i].B2c = Soln[j][i].B2c; pGrid->U[k][j][i].B3c = Soln[j][i].B3c; #ifndef ISOTHERMAL Soln[j][i].E = pres/Gamma_1 + 0.5*(SQR(pGrid->U[k][j][i].B1c) + SQR(pGrid->U[k][j][i].B2c) + SQR(pGrid->U[k][j][i].B3c) ) + 0.5*(SQR(pGrid->U[k][j][i].M1) + SQR(pGrid->U[k][j][i].M2) + SQR(pGrid->U[k][j][i].M3) )/den; pGrid->U[k][j][i].E = Soln[j][i].E; #endif } } } /* save solution on root grid */ if (pDomain->Level == 0) { for (j=js; j<=je; j++) { for (i=is; i<=ie; i++) { RootSoln[j][i].d = Soln[j][i].d ; RootSoln[j][i].M1 = Soln[j][i].M1; RootSoln[j][i].M2 = Soln[j][i].M2; RootSoln[j][i].M3 = Soln[j][i].M3; #ifndef ISOTHERMAL RootSoln[j][i].E = Soln[j][i].E ; #endif /* ISOTHERMAL */ #ifdef MHD RootSoln[j][i].B1c = Soln[j][i].B1c; RootSoln[j][i].B2c = Soln[j][i].B2c; RootSoln[j][i].B3c = Soln[j][i].B3c; #endif #if (NSCALARS > 0) for (n=0; n<NSCALARS; n++) RootSoln[j][i].s[n] = Soln[j][i].s[n]; #endif }} } return; }
void bvals_grav_init(MeshS *pM) { GridS *pG; DomainS *pD; int i,nl,nd,irefine; #ifdef MPI_PARALLEL int myL,myM,myN,l,m,n,nx1t,nx2t,nx3t,size; int x1cnt=0, x2cnt=0, x3cnt=0; /* Number of words passed in x1/x2/x3-dir. */ #endif /* MPI_PARALLEL */ /* Cycle through all the Domains that have active Grids on this proc */ for (nl=0; nl<(pM->NLevels); nl++){ for (nd=0; nd<(pM->DomainsPerLevel[nl]); nd++){ if (pM->Domain[nl][nd].Grid != NULL) { pD = (DomainS*)&(pM->Domain[nl][nd]); /* ptr to Domain */ pG = pM->Domain[nl][nd].Grid; /* ptr to Grid */ irefine = 1; for (i=1;i<=nl;i++) irefine *= 2; /* C pow fn only takes doubles !! */ #ifdef MPI_PARALLEL /* get (l,m,n) coordinates of Grid being updated on this processor */ get_myGridIndex(pD, myID_Comm_world, &myL, &myM, &myN); #endif /* MPI_PARALLEL */ /* Set function pointers for physical boundaries in x1-direction */ if(pG->Nx[0] > 1) { /*---- ix1 boundary ----------------------------------------------------------*/ if(pD->ix1_GBCFun == NULL){ /* Domain boundary is in interior of root */ if(pD->Disp[0] != 0) { pD->ix1_GBCFun = ProlongateLater; /* Domain is at L-edge of root Domain */ } else { switch(pM->BCFlag_ix1){ case 1: /* Reflecting */ pD->ix1_GBCFun = reflect_Phi_ix1; break; case 2: /* Outflow */ ath_error("[bvals_grav_init]: BCFlag_ix1 = 2 not implemented\n"); break; case 4: /* Periodic */ pD->ix1_GBCFun = periodic_Phi_ix1; #ifdef MPI_PARALLEL if(pG->lx1_Gid < 0 && pD->NGrid[0] > 1){ pG->lx1_Gid = pD->GData[myN][myM][pD->NGrid[0]-1].ID_Comm_Domain; } #endif /* MPI_PARALLEL */ break; case 5: /* Reflecting, B_normal!=0 */ pD->ix1_GBCFun = reflect_Phi_ix1; break; default: ath_error("[bvals_grav_init]: BCFlag_ix1 = %d unknown\n", pM->BCFlag_ix1); } } } /*---- ox1 boundary ----------------------------------------------------------*/ if(pD->ox1_GBCFun == NULL){ /* Domain boundary is in interior of root */ if((pD->Disp[0] + pD->Nx[0])/irefine != pM->Nx[0]) { pD->ox1_GBCFun = ProlongateLater; /* Domain is at R-edge of root Domain */ } else { switch(pM->BCFlag_ox1){ case 1: /* Reflecting */ pD->ox1_GBCFun = reflect_Phi_ox1; break; case 2: /* Outflow */ ath_error("[bvals_grav_init]: BCFlag_ox1 = 2 not implemented\n"); break; case 4: /* Periodic */ pD->ox1_GBCFun = periodic_Phi_ox1; #ifdef MPI_PARALLEL if(pG->rx1_Gid < 0 && pD->NGrid[0] > 1){ pG->rx1_Gid = pD->GData[myN][myM][0].ID_Comm_Domain; } #endif /* MPI_PARALLEL */ break; case 5: /* Reflecting, B_normal!=0 */ pD->ox1_GBCFun = reflect_Phi_ox1; break; default: ath_error("[bvals_grav_init]: BCFlag_ox1 = %d unknown\n", pM->BCFlag_ox1); } } } } /* Set function pointers for physical boundaries in x2-direction */ if(pG->Nx[1] > 1) { /*---- ix2 boundary ----------------------------------------------------------*/ if(pD->ix2_GBCFun == NULL){ /* Domain boundary is in interior of root */ if(pD->Disp[1] != 0) { pD->ix2_GBCFun = ProlongateLater; /* Domain is at L-edge of root Domain */ } else { switch(pM->BCFlag_ix2){ case 1: /* Reflecting */ pD->ix2_GBCFun = reflect_Phi_ix2; break; case 2: /* Outflow */ ath_error("[bvals_grav_init]: BCFlag_ix2 = 2 not implemented\n"); break; case 4: /* Periodic */ pD->ix2_GBCFun = periodic_Phi_ix2; #ifdef MPI_PARALLEL if(pG->lx2_Gid < 0 && pD->NGrid[1] > 1){ pG->lx2_Gid = pD->GData[myN][pD->NGrid[1]-1][myL].ID_Comm_Domain; } #endif /* MPI_PARALLEL */ break; case 5: /* Reflecting, B_normal!=0 */ pD->ix2_GBCFun = reflect_Phi_ix2; break; default: ath_error("[bvals_grav_init]: BCFlag_ix2 = %d unknown\n", pM->BCFlag_ix2); } } } /*---- ox2 boundary ----------------------------------------------------------*/ if(pD->ox2_GBCFun == NULL){ /* Domain boundary is in interior of root */ if((pD->Disp[1] + pD->Nx[1])/irefine != pM->Nx[1]) { pD->ox2_GBCFun = ProlongateLater; /* Domain is at R-edge of root Domain */ } else { switch(pM->BCFlag_ox2){ case 1: /* Reflecting */ pD->ox2_GBCFun = reflect_Phi_ox2; break; case 2: /* Outflow */ ath_error("[bvals_grav_init]: BCFlag_ox2 = 2 not implemented\n"); break; case 4: /* Periodic */ pD->ox2_GBCFun = periodic_Phi_ox2; #ifdef MPI_PARALLEL if(pG->rx2_Gid < 0 && pD->NGrid[1] > 1){ pG->rx2_Gid = pD->GData[myN][0][myL].ID_Comm_Domain; } #endif /* MPI_PARALLEL */ break; case 5: /* Reflecting, B_normal!=0 */ pD->ox2_GBCFun = reflect_Phi_ox2; break; default: ath_error("[bvals_grav_init]: BCFlag_ox2 = %d unknown\n", pM->BCFlag_ox2); } } } } /* Set function pointers for physical boundaries in x3-direction */ if(pG->Nx[2] > 1) { /*---- ix3 boundary ----------------------------------------------------------*/ if(pD->ix3_GBCFun == NULL){ /* Domain boundary is in interior of root */ if(pD->Disp[2] != 0) { pD->ix3_GBCFun = ProlongateLater; /* Domain is at L-edge of root Domain */ } else { switch(pM->BCFlag_ix3){ case 1: /* Reflecting */ pD->ix3_GBCFun = reflect_Phi_ix3; break; case 2: /* Outflow */ ath_error("[bvals_grav_init]: BCFlag_ix3 = 2 not defined\n"); break; case 4: /* Periodic */ pD->ix3_GBCFun = periodic_Phi_ix3; #ifdef MPI_PARALLEL if(pG->lx3_Gid < 0 && pD->NGrid[2] > 1){ pG->lx3_Gid = pD->GData[pD->NGrid[2]-1][myM][myL].ID_Comm_Domain; } #endif /* MPI_PARALLEL */ break; case 5: /* Reflecting, B_normal!=0 */ pD->ix3_GBCFun = reflect_Phi_ix3; break; default: ath_error("[bvals_grav_init]: BCFlag_ix3 = %d unknown\n", pM->BCFlag_ix3); } } } /*---- ox3 boundary ----------------------------------------------------------*/ if(pD->ox3_GBCFun == NULL){ /* Domain boundary is in interior of root */ if((pD->Disp[2] + pD->Nx[2])/irefine != pM->Nx[2]) { pD->ox3_GBCFun = ProlongateLater; /* Domain is at R-edge of root Domain */ } else { switch(pM->BCFlag_ox3){ case 1: /* Reflecting */ pD->ox3_GBCFun = reflect_Phi_ox3; break; case 2: /* Outflow */ ath_error("[bvals_grav_init]: BCFlag_ox3 = 2 not defined\n"); break; case 4: /* Periodic */ pD->ox3_GBCFun = periodic_Phi_ox3; #ifdef MPI_PARALLEL if(pG->rx3_Gid < 0 && pD->NGrid[2] > 1){ pG->rx3_Gid = pD->GData[0][myM][myL].ID_Comm_Domain; } #endif /* MPI_PARALLEL */ break; case 5: /* Reflecting, B_normal!=0 */ pD->ox3_GBCFun = reflect_Phi_ox3; break; default: ath_error("[bvals_grav_init]: BCFlag_ox3 = %d unknown\n", pM->BCFlag_ox3); } } } } /* Figure out largest size needed for send/receive buffers with MPI ----------*/ #ifdef MPI_PARALLEL for (n=0; n<(pD->NGrid[2]); n++){ for (m=0; m<(pD->NGrid[1]); m++){ for (l=0; l<(pD->NGrid[0]); l++){ /* x1cnt is surface area of x1 faces */ if(pD->NGrid[0] > 1){ nx2t = pD->GData[n][m][l].Nx[1]; if(nx2t > 1) nx2t += 1; nx3t = pD->GData[n][m][l].Nx[2]; if(nx3t > 1) nx3t += 1; if(nx2t*nx3t > x1cnt) x1cnt = nx2t*nx3t; } /* x2cnt is surface area of x2 faces */ if(pD->NGrid[1] > 1){ nx1t = pD->GData[n][m][l].Nx[0]; if(nx1t > 1) nx1t += 2*nghost; nx3t = pD->GData[n][m][l].Nx[2]; if(nx3t > 1) nx3t += 1; if(nx1t*nx3t > x2cnt) x2cnt = nx1t*nx3t; } /* x3cnt is surface area of x3 faces */ if(pD->NGrid[2] > 1){ nx1t = pD->GData[n][m][l].Nx[0]; if(nx1t > 1) nx1t += 2*nghost; nx2t = pD->GData[n][m][l].Nx[1]; if(nx2t > 1) nx2t += 2*nghost; if(nx1t*nx2t > x3cnt) x3cnt = nx1t*nx2t; } } }} #endif /* MPI_PARALLEL */ }}} /* End loop over all Domains with active Grids -----------------------*/ #ifdef MPI_PARALLEL /* Allocate memory for send/receive buffers and MPI_Requests */ size = x1cnt > x2cnt ? x1cnt : x2cnt; size = x3cnt > size ? x3cnt : size; size *= nghost; /* Multiply by the third dimension */ if (size > 0) { if((send_buf = (double**)calloc_2d_array(2,size,sizeof(double))) == NULL) ath_error("[bvals_init]: Failed to allocate send buffer\n"); if((recv_buf = (double**)calloc_2d_array(2,size,sizeof(double))) == NULL) ath_error("[bvals_init]: Failed to allocate recv buffer\n"); } if((recv_rq = (MPI_Request*) calloc_1d_array(2,sizeof(MPI_Request))) == NULL) ath_error("[bvals_init]: Failed to allocate recv MPI_Request array\n"); if((send_rq = (MPI_Request*) calloc_1d_array(2,sizeof(MPI_Request))) == NULL) ath_error("[bvals_init]: Failed to allocate send MPI_Request array\n"); #endif /* MPI_PARALLEL */ return; }
int main(int argc, char* argv[]) { /* argument variables */ int f1=0,f2=0,fi=1; char *defdir = "."; char *inbase = NULL, *postname = NULL; char *indir = defdir; /* fild variables */ FILE *fid; char fname[100]; /* data variables */ int i,*cpuid,*property,ntype; long p,n,*pid,*order; float header[20],**data,time[2],*typeinfo=NULL; /* Read Arguments */ for (i=1; i<argc; i++) { /* If argv[i] is a 2 character string of the form "-?" then: */ if(*argv[i] == '-' && *(argv[i]+1) != '\0' && *(argv[i]+2) == '\0') { switch(*(argv[i]+1)) { case 'd': /* -d <indir> */ indir = argv[++i]; break; case 'i': /* -i <basename> */ inbase = argv[++i]; break; case 's': /* -s <post-name> */ postname = argv[++i]; break; case 'f': /* -f <# range(f1:f2:fi)>*/ sscanf(argv[++i],"%d:%d:%d",&f1,&f2,&fi); if (f2 == 0) f2 = f1; break; case 'h': /* -h */ usage(argv[0]); break; default: usage(argv[0]); break; } } } /* Checkpoints */ if (inbase == NULL) sort_error("Please specify input file basename using -i option!\n"); if (postname == NULL) sort_error("Please specify posterior file name using -s option!\n"); if ((f1>f2) || (f2<0) || (fi<=0)) sort_error("Wrong number sequence in the -f option!\n"); /* ====================================================================== */ for (i=f1; i<=f2; i+=fi) { fprintf(stderr,"Processing file number %d...\n",i); /* Step 1: Read the file */ sprintf(fname,"%s/%s.%04d.%s.lis",indir,inbase,i,postname); fid = fopen(fname,"rb"); if (fid == NULL) sort_error("Fail to open output file %s!\n",fname); /* read header */ fread(header,sizeof(float),12,fid); fread(&ntype,sizeof(int),1,fid); if (i == f1) typeinfo = (float*)calloc_1d_array(ntype,sizeof(float)); fread(typeinfo,sizeof(float),ntype,fid); fread(&time,sizeof(float),2,fid); fread(&n,sizeof(long),1,fid); data = (float**)calloc_2d_array(n,7,sizeof(float)); property = (int*)calloc_1d_array(n,sizeof(int)); pid = (long*)calloc_1d_array(n,sizeof(long)); cpuid = (int*)calloc_1d_array(n,sizeof(int)); order = (long*)calloc_1d_array(n,sizeof(long)); /* read data */ for (p=0; p<n; p++) { fread(data[p],sizeof(float),7,fid); fread(&(property[p]),sizeof(int),1,fid); fread(&(pid[p]),sizeof(long),1,fid); fread(&(cpuid[p]),sizeof(int),1,fid); } fclose(fid); /* Step 2: sort the particles */ for (p=0; p<n; p++) order[p] = p; quicksort(cpuid, pid, order, 0, n-1); /* Step 3: write back the ordered particle list */ fid = fopen(fname,"wb"); /* write header */ fwrite(header,sizeof(float),12,fid); fwrite(&ntype,sizeof(int),1,fid); fwrite(typeinfo,sizeof(float),ntype,fid); fwrite(time,sizeof(float),2,fid); fwrite(&n,sizeof(long),1,fid); /* write data */ for (p=0; p<n; p++) { fwrite(data[order[p]],sizeof(float),7,fid); fwrite(&property[p],sizeof(int),1,fid); fwrite(&pid[order[p]],sizeof(long),1,fid); fwrite(&cpuid[order[p]],sizeof(int),1,fid); } free_2d_array(data); free(property); free(pid); free(cpuid); free(order); fclose(fid); } if (typeinfo != NULL) free(typeinfo); return 0; }
void problem(DomainS *pDomain) { GridS *pGrid = pDomain->Grid; int i,j,k,ks,pt,tsmode; long p,q; Real ScaleHg,tsmin,tsmax,tscrit,amin,amax,Hparmin,Hparmax; Real *ep,*ScaleHpar,epsum,mratio,pwind,rhoaconv,etavk; Real *epsilon,*uxNSH,*uyNSH,**wxNSH,**wyNSH; Real rhog,h,x1,x2,x3,t,x1p,x2p,x3p,zmin,zmax,dx3_1,b; long int iseed = myID_Comm_world; /* Initialize on the first call to ran2 */ if (pDomain->Nx[2] == 1) { ath_error("[par_strat3d]: par_strat3d only works for 3D problem.\n"); } #ifdef MPI_PARALLEL if (pDomain->NGrid[2] > 2) { ath_error( "[par_strat3d]: The z-domain can not be decomposed into more than 2 grids\n"); } #endif /* Initialize boxsize */ x1min = pGrid->MinX[0]; x1max = pGrid->MaxX[0]; Lx = x1max - x1min; x2min = pGrid->MinX[1]; x2max = pGrid->MaxX[1]; Ly = x2max - x2min; x3min = par_getd("domain1","x3min"); x3max = par_getd("domain1","x3max"); Lz = x3max - x3min; Lg = nghost*pGrid->dx3; /* size of the ghost zone */ ks = pGrid->ks; /* Read initial conditions */ Omega_0 = par_getd("problem","omega"); qshear = par_getd_def("problem","qshear",1.5); ipert = par_geti_def("problem","ipert",1); vsc1 = par_getd_def("problem","vsc1",0.05); /* in unit of iso_sound (N.B.!) */ vsc2 = par_getd_def("problem","vsc2",0.0); vsc1 = vsc1 * Iso_csound; vsc2 = vsc2 * Iso_csound; ScaleHg = Iso_csound/Omega_0; /* particle number */ Npar = (long)(par_geti("particle","parnumgrid")); pGrid->nparticle = Npar*npartypes; for (i=0; i<npartypes; i++) grproperty[i].num = Npar; if (pGrid->nparticle+2 > pGrid->arrsize) particle_realloc(pGrid, pGrid->nparticle+2); ep = (Real*)calloc_1d_array(npartypes, sizeof(Real)); ScaleHpar = (Real*)calloc_1d_array(npartypes, sizeof(Real)); epsilon = (Real*)calloc_1d_array(npartypes, sizeof(Real)); wxNSH = (Real**)calloc_2d_array(pGrid->Nx[2]+1, npartypes,sizeof(Real)); wyNSH = (Real**)calloc_2d_array(pGrid->Nx[2]+1, npartypes,sizeof(Real)); uxNSH = (Real*)calloc_1d_array(pGrid->Nx[2]+1, sizeof(Real)); uyNSH = (Real*)calloc_1d_array(pGrid->Nx[2]+1, sizeof(Real)); /* particle stopping time */ tsmode = par_geti("particle","tsmode"); if (tsmode == 3) {/* fixed stopping time */ tsmin = par_getd("problem","tsmin"); /* in code unit */ tsmax = par_getd("problem","tsmax"); tscrit= par_getd("problem","tscrit"); for (i=0; i<npartypes; i++) { tstop0[i] = tsmin*exp(i*log(tsmax/tsmin)/MAX(npartypes-1,1.0)); grproperty[i].rad = tstop0[i]; /* use fully implicit integrator for well coupled particles */ if (tstop0[i] < tscrit) grproperty[i].integrator = 3; } } else { amin = par_getd("problem","amin"); amax = par_getd("problem","amax"); for (i=0; i<npartypes; i++) grproperty[i].rad = amin*exp(i*log(amax/amin)/MAX(npartypes-1,1.0)); if (tsmode <= 2) {/* Epstein/General regime */ /* conversion factor for rhoa */ rhoaconv = par_getd_def("problem","rhoaconv",1.0); for (i=0; i<npartypes; i++) grrhoa[i]=grproperty[i].rad*rhoaconv; } if (tsmode == 1) /* General drag formula */ alamcoeff = par_getd("problem","alamcoeff"); } /* particle scale height */ Hparmax = par_getd("problem","hparmax"); /* in unit of gas scale height */ Hparmin = par_getd("problem","hparmin"); for (i=0; i<npartypes; i++) ScaleHpar[i] = Hparmax* exp(-i*log(Hparmax/Hparmin)/MAX(npartypes-1,1.0)); #ifdef FEEDBACK mratio = par_getd_def("problem","mratio",0.0); /* total mass fraction */ pwind = par_getd_def("problem","pwind",0.0); /* power law index */ if (mratio < 0.0) ath_error("[par_strat2d]: mratio must be positive!\n"); epsum = 0.0; for (i=0; i<npartypes; i++) { ep[i] = pow(grproperty[i].rad,pwind); epsum += ep[i]; } for (i=0; i<npartypes; i++) { ep[i] = mratio*ep[i]/epsum; grproperty[i].m = sqrt(2.0*PI)*ScaleHg/Lz*ep[i]* pGrid->Nx[0]*pGrid->Nx[1]*pGrid->Nx[2]/Npar; } #else mratio = 0.0; for (i=0; i<npartypes; i++) ep[i] = 0.0; #endif /* NSH equilibrium */ for (k=pGrid->ks; k<=pGrid->ke+1; k++) { h = pGrid->MinX[2] + (k-pGrid->ks)*pGrid->dx3; q = k - ks; etavk = fabs(vsc1+vsc2*SQR(h)); for (i=0; i<npartypes; i++) { epsilon[i] = ep[i]/ScaleHpar[i]*exp(-0.5*SQR(h/ScaleHg) *(SQR(1.0/ScaleHpar[i])-1.0))/erf(Lz/(sqrt(8.0)*ScaleHpar[i]*ScaleHg)); if (tsmode != 3) tstop0[i] = get_ts(pGrid,i,exp(-0.5*SQR(h/ScaleHg)),Iso_csound,etavk); } MultiNSH(npartypes, tstop0, epsilon, etavk, &uxNSH[q], &uyNSH[q], wxNSH[q], wyNSH[q]); } /* Now set initial conditions for the gas */ for (k=pGrid->ks; k<=pGrid->ke; k++) { for (j=pGrid->js; j<=pGrid->je; j++) { for (i=pGrid->is; i<=pGrid->ie; i++) { cc_pos(pGrid,i,j,k,&x1,&x2,&x3); rhog = exp(-0.5*SQR(x3/ScaleHg)); pGrid->U[k][j][i].d = rhog; if (ipert != 1) {/* NSH velocity */ pGrid->U[k][j][i].M1 = 0.5*rhog*(uxNSH[k-ks]+uxNSH[k-ks+1]); pGrid->U[k][j][i].M2 = 0.5*rhog*(uyNSH[k-ks]+uyNSH[k-ks+1]); } else { pGrid->U[k][j][i].M1 = 0.0; pGrid->U[k][j][i].M2 = 0.0; } pGrid->U[k][j][i].M3 = 0.0; #ifndef FARGO pGrid->U[k][j][i].M2 -= qshear*rhog*Omega_0*x1; #endif }}} /* Now set initial conditions for the particles */ p = 0; dx3_1 = 1.0/pGrid->dx3; zmin = pGrid->MinX[2]; zmax = pGrid->MaxX[2]; for (q=0; q<Npar; q++) { for (pt=0; pt<npartypes; pt++) { x1p = x1min + Lx*ran2(&iseed); x2p = x2min + Ly*ran2(&iseed); x3p = ScaleHpar[pt]*ScaleHg*Normal(&iseed); while ((x3p >= zmax) || (x3p < zmin)) x3p = ScaleHpar[pt]*ScaleHg*Normal(&iseed); pGrid->particle[p].property = pt; pGrid->particle[p].x1 = x1p; pGrid->particle[p].x2 = x2p; pGrid->particle[p].x3 = x3p; if (ipert != 1) {/* NSH velocity */ cellk(pGrid, x3p, dx3_1, &k, &b); k = k-pGrid->ks; b = b - pGrid->ks; pGrid->particle[p].v1 = (k+1-b)*wxNSH[k][pt]+(b-k)*wxNSH[k+1][pt]; pGrid->particle[p].v2 = (k+1-b)*wyNSH[k][pt]+(b-k)*wyNSH[k+1][pt]; } else { pGrid->particle[p].v1 = 0.0; pGrid->particle[p].v2 = vsc1+vsc2*SQR(x2p); } pGrid->particle[p].v3 = 0.0; #ifndef FARGO pGrid->particle[p].v2 -= qshear*Omega_0*x1p; #endif pGrid->particle[p].pos = 1; /* grid particle */ pGrid->particle[p].my_id = p; #ifdef MPI_PARALLEL pGrid->particle[p].init_id = myID_Comm_world; #endif p++; }} /* enroll gravitational potential function, shearing sheet BC functions */ ShearingBoxPot = StratifiedDisk; dump_history_enroll(hst_rho_Vx_dVy, "<rho Vx dVy>"); /* set the # of the particles in list output * (by default, output 1 particle per cell) */ nlis = par_geti_def("problem","nlis",pGrid->Nx[0]*pGrid->Nx[1]*pGrid->Nx[2]); /* set the number of particles to keep track of */ ntrack = par_geti_def("problem","ntrack",2000); /* set the threshold particle density */ dpar_thresh = par_geti_def("problem","dpar_thresh",10.0); /* finalize */ free(ep); free(ScaleHpar); free(epsilon); free_2d_array(wxNSH); free_2d_array(wyNSH); free(uxNSH); free(uyNSH); return; }