EXPORT int comps_consistent_at_node( NODE *node) { INTERFACE *intfc = node->interface; COMPONENT compi, compj; CURVE **c; int i, j; int ans; static O_NODE On; static int alloc_num_c = 0; /* Don't check for consistency on subdomain boundaries */ if (is_pp_node(node)) return YES; ans = YES; On.node = node; On.num_c = 0; On.prev = On.next = NULL; for (c = node->in_curves; c && *c; ++c) ++On.num_c; for (c = node->out_curves; c && *c; ++c) ++On.num_c; if (On.num_c == 0) return YES; if ((alloc_num_c == 0) || (alloc_num_c < On.num_c)) { if (alloc_num_c > 0) free_these(5,On.nc,On.nopp,On.pt,On.ang,On.orient); uni_array(&On.nc,On.num_c,sizeof(CURVE *)); uni_array(&On.nopp,On.num_c,sizeof(NODE *)); uni_array(&On.pt,On.num_c,sizeof(POINT *)); uni_array(&On.ang,On.num_c,FLOAT); uni_array(&On.orient,On.num_c,INT); alloc_num_c = On.num_c; } set_curves_at_onode(&On); for (i = 0; i < On.num_c; ++i) { j = (i + 1) % On.num_c; if (On.orient[i] == POSITIVE_ORIENTATION) compi = negative_component(On.nc[i]); else compi = positive_component(On.nc[i]); if (On.orient[j] == POSITIVE_ORIENTATION) compj = positive_component(On.nc[j]); else compj = negative_component(On.nc[j]); if (!equivalent_comps(compi,compj,intfc)) { ans = NO; break; } } return ans; } /*end comps_consistent_at_node*/
EXPORT const COMPONENT *comps_with_type( COMP_TYPE_TYPE type, int *num) { static COMPONENT *compstore = NULL; static int len_compstore = 0; int i, N; if (len_compstore <= max_n_comps) { if (compstore != NULL) free(compstore); len_compstore = max_n_comps+1; uni_array(&compstore,len_compstore,sizeof(COMPONENT)); } for (i = 0; i <= max_n_comps; ++i) compstore[i] = NO_COMP; for (N = 0, i = 0; i <= max_n_comps; ++i) { if ((ct[i] != NULL) && (ct[i]->type == type)) compstore[N++] = ct[i]->comp; } if (num != NULL) *num = N; return (N > 0) ? compstore : NULL; } /*end comps_with_type*/
EXPORT COMP_TYPE *comp_type( COMPONENT comp) { static int maxcomp = -1; if (ct == NULL) { int i; uni_array(&ct,max_n_comps+1,sizeof(COMP_TYPE*)); for (i = 0; i <= max_n_comps; ++i) ct[i] = NULL; } if (comp < 0) comp = maxcomp + 1; if (comp > max_n_comps) { COMP_TYPE **new_ct; int i, new_max_n_comps = max(comp,2*max_n_comps); uni_array(&new_ct,new_max_n_comps+1,sizeof(COMP_TYPE*)); for (i = 0; i <= max_n_comps; ++i) new_ct[i] = ct[i]; for (;i <= new_max_n_comps; ++i) new_ct[i] = NULL; free(ct); ct = new_ct; max_n_comps = new_max_n_comps; } if (maxcomp < comp) maxcomp = comp; if (ct[comp] == NULL) ct[comp] = new_comp_type(comp); else if (ct[comp]->comp != comp) { screen("\nERROR in comp_type(), Inconsistent component number!\n"); clean_up(ERROR); } return ct[comp]; } /*end comp_type*/
EXPORT int set_grid_lines( RECT_GRID *rgr) { register int j, jmin, jmax; register double *edges, *centers, *dh, h, L; double *store; int i, dim = rgr->dim; int lbuf[MAXD], ubuf[MAXD]; int len[MAXD], tlen; static const int EXTRA_BUF = 2; if (rgr == NULL) { return NO; } tlen = 0; for (i = 0; i < dim; ++i) { lbuf[i] = rgr->lbuf[i] + EXTRA_BUF; ubuf[i] = rgr->ubuf[i] + EXTRA_BUF; len[i] = (lbuf[i] + rgr->gmax[i] + ubuf[i]); tlen += 3*len[i]+1; } uni_array(&rgr->glstore,tlen,FLOAT); if (rgr->glstore == NULL) { return NO; } store = rgr->glstore; for (i = 0; i < dim; ++i) { rgr->edges[i] = store + lbuf[i]; store += len[i]+1; rgr->centers[i] = store + lbuf[i]; store += len[i]; rgr->dh[i] = store + lbuf[i]; store += len[i]; jmin = -lbuf[i]; jmax = rgr->gmax[i] + ubuf[i]; L = rgr->L[i]; h = rgr->h[i]; edges = rgr->edges[i] + jmin; centers = rgr->centers[i] + jmin; dh = rgr->dh[i] + jmin; for (j = jmin; j < jmax; ++j) { *edges = L + j*h; *centers++ = *edges++ + 0.5*h; *dh++ = h; } *edges = L + jmax*h; } return YES; } /*end set_grid_lines*/
LOCAL char *set_ppfname( char *ppfname, const char *fname, size_t *ppfname_len) { size_t len; if (pp_numnodes() > 1 ) { const char *nd; nd = right_flush(pp_mynode(),4); if (fname == NULL) fname = ""; len = strlen(fname)+1+strlen(nd)+1; if (*ppfname_len < len) { *ppfname_len = len; if (ppfname != NULL) free(ppfname); uni_array(&ppfname,*ppfname_len,CHAR); } (void) sprintf(ppfname,"%s.%s",fname,nd); } else { if (fname == NULL) fname = ""; len = strlen(fname)+1; if (*ppfname_len < len) { *ppfname_len = len; if (ppfname != NULL) free(ppfname); uni_array(&ppfname,*ppfname_len,CHAR); } (void) strcpy(ppfname,fname); } return ppfname; } /*end set_ppfname*/
EXPORT bool FrontIntrpStateAtPoint( Front *front, COMPONENT comp, float *coords, float *grid_array, float (*get_state)(Locstate), float *ans) { int icoords[MAXD]; INTERFACE *grid_intfc = front->grid_intfc; static INTRP_CELL *blk_cell; RECT_GRID *gr = &topological_grid(grid_intfc); int i,dim = gr->dim; if (blk_cell == NULL) { scalar(&blk_cell,sizeof(INTRP_CELL)); uni_array(&blk_cell->var,MAX_NUM_VERTEX_IN_CELL,sizeof(float)); bi_array(&blk_cell->coords,MAX_NUM_VERTEX_IN_CELL,MAXD,sizeof(float)); bi_array(&blk_cell->p_lin,MAXD+1,MAXD,sizeof(float)); uni_array(&blk_cell->var_lin,MAXD+1,sizeof(float)); lin_cell_tol = 1.0; for (i = 0; i < dim; ++i) lin_cell_tol *= 0.00001*gr->h[i]; } if (!rect_in_which(coords,icoords,gr)) { return NO; } collect_cell_ptst(blk_cell,icoords,comp,front,grid_array,get_state); if (build_linear_element(blk_cell,coords)) { *ans = FrontLinIntrp(coords,blk_cell,NO); return YES; } else return NO; } /* end FrontIntrpStateAtPoint */
LOCAL char *get_list_file_name( char *fname, const char *dname, const char *name, size_t *fname_len) { size_t len; if (dname == NULL) dname = ""; len = strlen(dname)+1+strlen(name)+6; if (*fname_len < len) { *fname_len = len; if (fname != NULL) free(fname); uni_array(&fname,*fname_len,CHAR); } if (strlen(dname) != 0) (void) sprintf(fname,"%s/%s.list",dname,name); else (void) sprintf(fname,"%s.list",name); return fname; }
EXPORT PP_GRID *set_pp_grid( INIT_DATA *init, RECT_GRID *comp_glbgr) { float L[MAXD], U[MAXD], *GL, *GU; float *h = comp_glbgr->h; int lbuf[MAXD], ubuf[MAXD]; int gmax[MAXD]; int icoords[MAXD]; int i, dim = comp_glbgr->dim; int myid = pp_mynode(); static PP_GRID Pp_grid; PP_GRID *pp_grid = &Pp_grid; char s[1028]; debug_print("init_pp_grid","Entered set_pp_grid():\n"); copy_rect_grid(&pp_grid->Global_grid,comp_glbgr); pp_grid->nn = 1; for (i = 0; i < dim; ++i) { int Gmax, Pmax, k; int basic_slices, extra_slices; pp_grid->buf[i] = buffer_zones(init)[i]; Pmax = pp_grid->gmax[i] = subdomains(init)[i]; pp_grid->nn *= Pmax; uni_array(&pp_grid->dom[i],Pmax + 1,FLOAT); pp_grid->dom[i][0] = comp_glbgr->L[i]; pp_grid->dom[i][Pmax] = comp_glbgr->U[i]; Gmax = comp_glbgr->gmax[i]; basic_slices = Gmax / Pmax; extra_slices = Gmax % Pmax; for (k = 1; k < Pmax; ++k) { if (k < extra_slices) pp_grid->dom[i][k] = k*(basic_slices + 1)*h[i] + pp_grid->dom[i][0]; else pp_grid->dom[i][k] = (k*basic_slices + extra_slices)*h[i] + pp_grid->dom[i][0]; } } /* Clip rectangular grid to subdomain */ GL = pp_grid->Global_grid.L; GU = pp_grid->Global_grid.U; find_Cartesian_coordinates(myid,pp_grid,icoords); for (i = 0; i < dim; ++i) { L[i] = pp_grid->dom[i][icoords[i]]; U[i] = pp_grid->dom[i][icoords[i] + 1]; gmax[i] = irint((U[i] - L[i])/h[i]); switch (dim) /* TODO Unify 2 and 3 D */ { case 1: case 2: lbuf[i] = (icoords[i] > 0) ? pp_grid->buf[i] : 0; ubuf[i] = (icoords[i]<(pp_grid->gmax[i]-1))?pp_grid->buf[i]:0; break; case 3: lbuf[i] = pp_grid->buf[i]; ubuf[i] = pp_grid->buf[i]; break; } } set_rect_grid(L,U,GL,GU,lbuf,ubuf,gmax,dim,&comp_glbgr->Remap, &pp_grid->Zoom_grid); if (debugging("init_pp_grid")) { (void) printf("pp_grid after set_pp_grid()\n"); (void) print_PP_GRID_structure(pp_grid); } debug_print("init_pp_grid","Left i_set_pp_grid():\n"); return pp_grid; } /*end set_pp_grid*/
EXPORT bool output_spectral_analysis( char *basename, Wave *wave, Front *front) { COMPONENT comp; Locstate state; float *coords; float *L = wave->rect_grid->L; float *U = wave->rect_grid->U; float *h = wave->rect_grid->h; int icoords[MAXD]; int dim = wave->rect_grid->dim; int status; int step = front->step; char energy_name[100],vorticity_name[100], enstrophy_name[100],dens_name[100],pres_name[100]; FILE *energy_file,*vorticity_file, *enstrophy_file,*dens_file,*pres_file; debug_print("fft","Entered fft_energy_spectral()\n"); (void) sprintf(energy_name,"%s.energy%s.dat",basename, right_flush(step,TSTEP_FIELD_WIDTH)); (void) sprintf(vorticity_name,"%s.vorticity%s.dat",basename, right_flush(step,TSTEP_FIELD_WIDTH)); (void) sprintf(enstrophy_name,"%s.enstrophy%s.dat",basename, right_flush(step,TSTEP_FIELD_WIDTH)); (void) sprintf(dens_name,"%s.density%s.dat",basename, right_flush(step,TSTEP_FIELD_WIDTH)); (void) sprintf(pres_name,"%s.pressure%s.dat",basename, right_flush(step,TSTEP_FIELD_WIDTH)); energy_file = fopen(energy_name,"w"); vorticity_file = fopen(vorticity_name,"w"); enstrophy_file = fopen(enstrophy_name,"w"); dens_file = fopen(dens_name,"w"); pres_file = fopen(pres_name,"w"); if (wave->sizest == 0) { debug_print("fft","Left fft_energy_spectral()\n"); return FUNCTION_FAILED; } switch (dim) { #if defined(ONED) case 1: { int ix; int xmax; COMPLEX *mesh_energy; xmax = wave->rect_grid->gmax[0]; uni_array(&mesh_energy,xmax,sizeof(COMPLEX)); for (ix = 0; ix < xmax; ++ix) { icoords[0] = ix; coords = Rect_coords(icoords,wave); comp = Rect_comp(icoords,wave); state = Rect_state(icoords,wave); mesh_energy[ix].real = Energy(state); mesh_energy[ix].imag = 0.0; } break; } #endif /* defined(ONED) */ #if defined(TWOD) case 2: { int ix, iy; int xmax, ymax, mx, my, dummy; COMPLEX **mesh_energy,**mesh_vorticity,**mesh_enstrophy; Locstate lstate,rstate,bstate,tstate; float kk,kx,ky,dk; xmax = wave->rect_grid->gmax[0]; ymax = wave->rect_grid->gmax[1]; if (!Powerof2(xmax,&mx,&dummy) || !Powerof2(ymax,&my,&dummy)) { screen("fft_energy_spectral() cannot analyze " "mesh not power of 2\n"); screen("xmax = %d ymax = %d\n",xmax,ymax); return FUNCTION_FAILED; } bi_array(&mesh_energy,xmax,ymax,sizeof(COMPLEX)); bi_array(&mesh_vorticity,xmax,ymax,sizeof(COMPLEX)); fprintf(energy_file,"zone i=%d, j=%d\n",xmax,ymax); fprintf(vorticity_file,"zone i=%d, j=%d\n",xmax,ymax); fprintf(dens_file,"zone i=%d, j=%d\n",xmax,ymax); fprintf(pres_file,"zone i=%d, j=%d\n",xmax,ymax); fprintf(enstrophy_file,"zone i=%d, j=%d\n",xmax,ymax); for (iy = 0; iy < ymax; ++iy) { for (ix = 0; ix < xmax; ++ix) { icoords[0] = ix; icoords[1] = iy; coords = Rect_coords(icoords,wave); comp = Rect_comp(icoords,wave); state = Rect_state(icoords,wave); mesh_energy[ix][iy].real = kinetic_energy(state); mesh_energy[ix][iy].imag = 0.0; if (ix != 0) icoords[0] = ix - 1; else icoords[0] = ix; lstate = Rect_state(icoords,wave); if (ix != xmax-1) icoords[0] = ix + 1; else icoords[0] = ix; rstate = Rect_state(icoords,wave); icoords[0] = ix; if (iy != 0) icoords[1] = iy - 1; else icoords[1] = iy; bstate = Rect_state(icoords,wave); if (iy != ymax-1) icoords[1] = iy + 1; else icoords[1] = iy; tstate = Rect_state(icoords,wave); mesh_vorticity[ix][iy].real = (Mom(rstate)[1]/Dens(rstate) - Mom(lstate)[1]/Dens(lstate))/h[0] - (Mom(tstate)[0]/Dens(tstate) - Mom(bstate)[0]/Dens(bstate))/h[1]; mesh_vorticity[ix][iy].imag = 0.0; fprintf(energy_file,"%lf\n",kinetic_energy(state)); fprintf(vorticity_file,"%lf\n",mesh_vorticity[ix][iy].real); fprintf(enstrophy_file,"%lf\n", sqr(mesh_vorticity[ix][iy].real)); fprintf(dens_file,"%lf\n",Dens(state)); fprintf(pres_file,"%lf\n",pressure(state)); } } fft_output2d(basename,"energy",step,wave->rect_grid, mesh_energy); fft_output2d(basename,"vorticity",step,wave->rect_grid, mesh_vorticity); free_these(2,mesh_energy,mesh_vorticity); break; } #endif /* defined(TWOD) */ #if defined(THREED) case 3: { int ix, iy, iz; int xmax, ymax, zmax; COMPLEX ***mesh_energy; xmax = wave->rect_grid->gmax[0]; ymax = wave->rect_grid->gmax[1]; zmax = wave->rect_grid->gmax[2]; tri_array(&mesh_energy,xmax,ymax,zmax,sizeof(COMPLEX)); for (iz = 0; iz < zmax; ++iz) { icoords[2] = iz; for (iy = 0; iy < ymax; ++iy) { icoords[1] = iy; for (ix = 0; ix < xmax; ++ix) { icoords[0] = ix; coords = Rect_coords(icoords,wave); comp = Rect_comp(icoords,wave); state = Rect_state(icoords,wave); mesh_energy[ix][iy][iz].real = Energy(state); mesh_energy[ix][iy][iz].imag = 0.0; } } } break; } #endif /* defined(THREED) */ } fclose(energy_file); fclose(vorticity_file); fclose(enstrophy_file); fclose(dens_file); fclose(pres_file); debug_print("fft","Left fft_energy_spectral()\n"); return FUNCTION_SUCCEEDED; } /*end fft_energy_spectral*/
/*ARGSUSED*/ EXPORT void u_pp_send( int tag, POINTER buf, size_t len, int node, const char *file, int line) { #if defined(__MPI__) static byte *msg_buf = NULL; int mpi_return_status; #endif /* defined(__MPI__) */ if (debugging("pp_clock")) start_clock("pp_send"); if (DEBUG) { (void) printf("Node %d sending message with tag %d ", pp_mynode(),tag); (void) printf("and len %d to node %d, ",(int)len,node); (void) printf("File %s, line %d\n",file,line); } pp_okay_to_proceed("u_pp_send","no message sent"); #if defined(__MPI__) if (msg_buf == NULL) { uni_array(&msg_buf,MSG_BUF_SIZE,sizeof(byte)); mpi_return_status = MPI_Buffer_attach(msg_buf,(int)MSG_BUF_SIZE); (void) printf("In first call to u_pp_send(), "); (void) printf("setting the buffer size to %lu bytes.\n", MSG_BUF_SIZE); if (mpi_return_status != MPI_SUCCESS) { screen("ERROR in u_pp_send(), " "MPI_Buffer_attach failed, " "mpi_return_status = %d\n",mpi_return_status); clean_up(ERROR); } } mpi_return_status = MPI_Bsend(buf,(int)len,MPI_BYTE, node,tag,FronTier_COMM); if (mpi_return_status != MPI_SUCCESS) { screen("ERROR in u_pp_send(), MPI_Send() failed, " "mpi_return_status = %d\n",mpi_return_status); clean_up(ERROR); } #else /* defined(__MPI__) */ COMMAND_NOT_IMPLEMENTED("u_pp_send","scalar mode"); #endif /* defined(__MPI__) */ if (debugging("pp_clock")) stop_clock("pp_send"); DEBUG_LEAVE(u_pp_send) } /*end u_pp_send*/
EXPORT bool make_point_comp_lists( INTERFACE *intfc) { int ix, ix0, ix1; int zp; POINT **p; struct Table *T; COMPONENT icomp; RECT_GRID *grid; if (DEBUG) (void) printf("Entered make_point_comp_lists()\n"); if ((T = table_of_interface(intfc)) == NULL) { (void) printf("WARNING in make_point_comp_lists(), " "table_of_interface = NULL\n"); return FUNCTION_FAILED; } if (no_topology_lists(intfc) == YES) { screen("ERROR in make_point_comp_lists(), " "illegal attempt to construct interface topology\n" "no_topology_lists(intfc) == YES\n"); clean_up(ERROR); } grid = &T->rect_grid; /* Free old storage if necessary */ if (T->num_of_points != NULL) free(T->num_of_points); if (T->pts_in_zone != NULL) free(T->pts_in_zone); if (T->compon1d != NULL) free(T->compon1d); /* Create a Grid if Needed: */ if (!T->fixed_grid) set_topological_grid(intfc,(RECT_GRID *)NULL); /* Allocate new storage if necessary */ uni_array(&T->num_of_points,grid->gmax[0],INT); if (T->num_of_points == NULL) { (void) printf("WARNING in make_point_complist(), " "can't allocate T->num_of_points\n"); return FUNCTION_FAILED; } if (DEBUG) (void) printf("T->num_of_points allocated\n"); /* NOTE: vector returns integer values initalized to zero */ uni_array(&T->compon1d,grid->gmax[0],sizeof(COMPONENT)); if (T->compon1d == NULL) { (void) printf("WARNING in make_point_complist(), " "can't allocate T->compon1d\n"); return FUNCTION_FAILED; } if (DEBUG) (void) printf("T->compon1d allocated\n"); for (ix = 0; ix < grid->gmax[0]; ++ix) T->compon1d[ix] = NO_COMP; uni_array(&T->pts_in_zone,grid->gmax[0],sizeof(POINT **)); /* NOTE: vector returns pointer values initalized to NULL */ if (T->pts_in_zone == NULL) { (void) printf("WARNING in make_point_complist(), " "can't allocate T->pts_in_zone\n"); return FUNCTION_FAILED; } start_clock("make_point_comp_lists"); if ((intfc->modified) && (intfc->num_points > 0)) sort_point_list(intfc->points,intfc->num_points); /* Default setting */ ix0 = 0; if (intfc->num_points != 0) { p = intfc->points; icomp = negative_component(*p); ix1 = -1; for (p = intfc->points; p && *p; ++p) { if (rect_in_which(Coords(*p),&zp,&T->rect_grid) == FUNCTION_FAILED) continue; ++T->num_of_points[zp]; if ((p == intfc->points) || (zp != ix1)) { T->pts_in_zone[zp] = p; } T->compon1d[zp] = ONFRONT; ix1 = zp; for (ix = ix0; ix < ix1; ++ix) T->compon1d[ix] = icomp; icomp = positive_component(*p); ix0 = zp + 1; } for (ix = ix0; ix < grid->gmax[0]; ++ix) T->compon1d[ix] = icomp; } stop_clock("make_point_comp_lists"); intfc->modified = NO; intfc->table->new_grid = NO; if (DEBUG) (void) printf("Leaving make_point_comp_lists()\n\n"); return FUNCTION_SUCCEEDED; } /*end make_point_comp_lists*/
//FIXME: Put me somewhere double blasius_f( double eta) { // return 0.5; // There exists a preexisting table relating eta values to f values // for the blasius function. They will be stored in these arrays static double *table_eta = NULL; static double *table_f = NULL; static int table_size; // Read in tabular blasius function // format is: // num_values(int) // eta(double) f(double) // ... ... if(table_eta == NULL) { FILE *f = fopen("blas.data","rt"); if(!f) { screen("Cannot find blas.data to initialize boundary" "layer.\nClosing cleanly.\n"); clean_up(-1); } int i, numi; fscanf(f,"%d\n",&numi); table_size = numi; uni_array(&table_eta, numi, sizeof(double)); uni_array(&table_f, numi, sizeof(double)); for(i=0; i < numi; ++i) { // fscanf(f, "%f %f\n",&table_eta[i], &table_f[i]); char buf[32]; fscanf(f,"%s",buf); table_eta[i] = atof(buf); fscanf(f,"%s",buf); table_f[i] = atof(buf); printf("%f %f\n",table_eta[i], table_f[i]); } fclose(f); } // Do bisection method to find which tabular indices // the desired point is on. int a=0, b=table_size-1; int g = (a+b)/2; while(b-a>1) { if(table_eta[g] < eta) { a=g; } else if(table_eta[g] > eta) { b=g; } else { return table_eta[g]; } g = (a+b)/2; // FIXME: prevent infinite loop in case of error } double f1 = table_f[a], f2 = table_f[b]; double e1 = table_eta[a], e2 = table_eta[b]; // Do linear interpolation from the table data to the known point. double coeff = (eta-e1)/(e2-e1); // printf("coeff=%f\n",coeff); double f= coeff*(f1) + (1.-coeff)*f2; // printf("blas_f=%f\n",f); return f; }
EXPORT void geomview_amr_fronts_plot2d( const char *dname, Front **frs, int num_patches, Wv_on_pc **redistr_table, int max_n_patch) { FILE *fp; int i; char fmt[256]; static const char *indent = " "; static char *fname = NULL, *ppfname = NULL; static size_t fname_len = 0, ppfname_len = 0; INTERFACE *intfc = frs[0]->interf; INTERFACE *tmpintfc; CURVE **c; INTERFACE *sav_intfc; bool sav_copy; float **clrmap = NULL; float ccolor[4] = {0.0, 0.0, 0.0, 1.0}; int myid, numnodes; const char *nstep; char outname[256],outdir[256]; myid = pp_mynode(); numnodes = pp_numnodes(); sprintf(outdir,"%s/%s",dname,"geomintfc"); ppfname = set_ppfname(ppfname,"intfc",&ppfname_len); nstep = right_flush(frs[0]->step,7); sprintf(outname,"%s.ts%s",ppfname,nstep); if (create_directory(dname,YES) == FUNCTION_FAILED) { (void) printf("WARNING in geomview_intfc_plot2d(), directory " "%s doesn't exist and can't be created\n",dname); return; } if (create_directory(outdir,YES) == FUNCTION_FAILED) { (void) printf("WARNING in geomview_intfc_plot2d(), directory " "%s doesn't exist and can't be created\n",outdir); return; } sav_intfc = current_interface(); sav_copy = copy_intfc_states(); set_size_of_intfc_state(size_of_state(intfc)); set_copy_intfc_states(YES); tmpintfc = copy_interface(intfc); /* clip_to_interior_region(tmpintfc, frs[0]->rect_grid->lbuf,frs[0]->rect_grid->ubuf); */ uni_array(&clrmap,6,sizeof(float*)); for(i = 0; i < 6; i++) uni_array(&clrmap[i],4,sizeof(float)); i = 0; clrmap[i][0] = 0.098; clrmap[i][1] = 0.647; clrmap[i][2] = 0.400; clrmap[i][3] = 1.000; i++; clrmap[i][0] = 0.898; clrmap[i][1] = 0.400; clrmap[i][2] = 0.000; clrmap[i][3] = 1.000; i++; clrmap[i][0] = 0.500; clrmap[i][1] = 1.000; clrmap[i][2] = 0.500; clrmap[i][3] = 1.000; i++; clrmap[i][0] = 1.000; clrmap[i][1] = 0.000; clrmap[i][2] = 1.000; clrmap[i][3] = 1.000; i++; clrmap[i][0] = 0.000; clrmap[i][1] = 0.800; clrmap[i][2] = 1.000; clrmap[i][3] = 1.000; i++; clrmap[i][0] = 0.250; clrmap[i][1] = 0.808; clrmap[i][2] = 0.098; clrmap[i][3] = 1.000; i++; fname = get_list_file_name(fname,outdir,outname,&fname_len); if ((fp = fopen(fname,"w")) == NULL) { (void) printf("WARNING in gview_plot_intfc2d(), " "can't open %s\n",fname); delete_interface(tmpintfc); set_current_interface(sav_intfc); set_copy_intfc_states(sav_copy); return; } fprintf(fp,"{ LIST \n"); /* beginning of writting Vect to file */ for(c = tmpintfc->curves; c and *c; c++) gview_plot_curve2d(fp,*c,ccolor); for(int i = 0; i < num_patches; i++) { int use_clr; /* gview_plot_box2d(fp, frs[i]->rect_grid->L, frs[i]->rect_grid->U,clrmap[i%3]); */ if(NULL != redistr_table) use_clr = redistr_table[myid][i].pc_id % 6; else use_clr = 1; gview_plot_grid2d(fp,frs[i]->rect_grid,clrmap[use_clr]); } /* end of LIST OBJ */ fprintf(fp,"}\n"); fclose(fp); set_current_interface(sav_intfc); set_copy_intfc_states(sav_copy); delete_interface(tmpintfc); for(int i = 0; i < 6; i++) free(clrmap[i]); free(clrmap); }
EXPORT int collect_pcs_in_mesh3d(TRI_GRID *ntg) { struct Table *T; P_LINK *hash_table; Locstate *states; COMPONENT *comp; TG_PT *node_pts; BLK_EL0 *blk_el0; TRI **tris; SURFACE **surfs; POINT_COMP_ST *pcs; int ix, iy, iz, l, nt, ***num_tris;; int *offset = ntg->offset; int num_pcs = 0; int h_size; int xmax, ymax, zmax; xmax = ntg->rect_grid.gmax[0]; ymax = ntg->rect_grid.gmax[1]; zmax = ntg->rect_grid.gmax[2]; ntg->_locate_on_trigrid = tg_build; set_tri3d_tolerances(ntg); T = table_of_interface(ntg->grid_intfc); //orig_construct_tri_grid //reconstruct_intfc_and_tri_grid // init_triangulation_storage // components, states // set_crx_structure_storage3dv0 // ntg->n_bilin_els = xmax*ymax*zmax; //expanded_topological_grid // ntg->n_node_points = (xmax+1)*(ymax+1)*(zmax+1); // // alloc_node_points(ntg,ntg->n_node_points); // node_points // alloc_blk_els0(ntg,ntg->n_bilin_els); // blk_els0 // set_interpolation_storage3dv0 // count_num_pcs3d // n_pcs, pcs, front_points /* Allocate space for hashing table */ h_size = (ntg->grid_intfc->num_points)*4+1; uni_array(&hash_table,h_size,sizeof(P_LINK)); copy_tg_pts_from_intfc(ntg,hash_table,h_size); //front_points copy_tg_pts_from_regular_grid(ntg); //node_points /* Triangulate each mesh block */ comp = T->components; states = ntg->states; node_pts = ntg->node_points; blk_el0 = ntg->blk_els0; num_tris = T->num_of_tris; for (iz = 0; iz < zmax; ++iz) { for (iy = 0; iy < ymax; ++iy) { for (ix = 0; ix < xmax; ++ix) { //debug with tg_build line_pj remove_from_debug("pcs_cell"); //if((ix == 3 && iy == 5 && iz == 0 && pp_mynode() == 0) || // (ix == 33 && iy == 5 && iz == 0 && pp_mynode() == 2) ) // add_to_debug("pcs_cell"); nt = num_tris[iz][iy][ix]; pcs = blk_el0_pcs_els(blk_el0) = &(ntg->pcs[num_pcs]); if (nt != 0) { tris = T->tris[iz][iy][ix]; surfs = T->surfaces[iz][iy][ix]; if(debugging("pcs_cell")) { int nd; printf("#pcs_cell %d\n", nt); for(nd = 0; nd < nt; nd++) { print_tri(tris[nd], surfs[nd]->interface); printf("%d (%d %d)\n", surfs[nd], negative_component(surfs[nd]), positive_component(surfs[nd])); } } for (l = 0; l < 8; ++l) { pcs[l].p = node_pts + offset[l]; pcs[l].comp[0] = comp[offset[l]]; pcs[l].s[0] = states[offset[l]]; pcs[l].comp[1] = NO_COMP; pcs[l].s[1] = NULL; } num_pcs_els_in_blk(blk_el0) = 8; add_intfc_blk_pcs3d(nt,tris,surfs,blk_el0, hash_table,h_size); num_pcs += num_pcs_els_in_blk(blk_el0); } else { set_bilinear_blk_el0(blk_el0); pcs[0].p = node_pts; pcs[0].comp[0] = comp[0]; pcs[0].s[0] = states[0]; pcs[0].comp[1] = NO_COMP; pcs[0].s[1] = NULL; ++num_pcs; } ++node_pts; ++states; ++comp; ++blk_el0; } ++node_pts; ++states; ++comp; } if (iz < zmax-1) { //node_pts, states, comp are defined on the nodes of the topological grid. node_pts += xmax+1; states += xmax+1; comp += xmax+1; } } //printf("#collect_pcs_in_mesh3d af: num_pcs = %d ntg->n_pcs=%d\n", // num_pcs, ntg->n_pcs); if(num_pcs != ntg->n_pcs) { printf("ERROR: num_pcs != ntg->n_pcs\n"); clean_up(ERROR); } free(hash_table); return GOOD_STEP; } /*end collect_pcs_in_mesh3d*/
EXPORT int collect_pcs_in_mesh2d(TRI_GRID *ntg) { int xmax; int ymax; struct Table *T; P_LINK *hash_table; Locstate *states; COMPONENT *comp; TG_PT *node_pts; BLK_EL0 *blk_el0; BOND ****by; BOND ***byx; CURVE ****cy,***cyx; POINT_COMP_ST *pcs; int ix, iy, l; int **nby, *nbyx; int *offset = ntg->offset; int num_pcs = 0; int h_size; #if defined(DEBUG_TRI_GRID) debug_print("collect_pcs_in_mesh2d", "make_bond_comp_list"); #endif /* defined(DEBUG_TRI_GRID) */ set_dual_interface_topology(ntg); ntg->_locate_on_trigrid = tg_build; #if defined(DEBUG_TRI_GRID) debug_print("collect_pcs_in_mesh2d", "to count_num_pcs2d"); #endif /* defined(DEBUG_TRI_GRID) */ ntg->n_pcs = count_num_pcs2d(ntg); VECTOR(ntg,pcs,ntg->n_pcs,sizeof(POINT_COMP_ST)); T = table_of_interface(ntg->grid_intfc); xmax = ntg->rect_grid.gmax[0]; ymax = ntg->rect_grid.gmax[1]; /* Allocate space for hashing table */ /*P_LINK (int.h) pair of left and right states */ h_size = (ntg->grid_intfc->num_points)*4+1; uni_array(&hash_table,h_size,sizeof(P_LINK)); #if defined(DEBUG_TRI_GRID) debug_print("collect_pcs_in_mesh2d", "to copy_tg_pts_from_intfc"); #endif /* defined(DEBUG_TRI_GRID) */ copy_tg_pts_from_intfc(ntg,hash_table,h_size); /* set the pcs's */ comp = T->components; states = ntg->states; node_pts = ntg->node_points; blk_el0 = ntg->blk_els0; #if defined(DEBUG_TRI_GRID) debug_print("collect_pcs_in_mesh2d", "to set the pcs's"); #endif /* defined(DEBUG_TRI_GRID) */ for (iy = 0, nby = T->num_of_bonds, by = T->bonds, cy = T->curves; iy < ymax; ++iy, ++nby, ++by, ++cy) { for (ix = 0, nbyx = *nby, byx = *by, cyx = *cy; ix < xmax; ++ix, ++nbyx, ++byx, ++cyx) { pcs = blk_el0_pcs_els(blk_el0) = &(ntg->pcs[num_pcs]); if (*nbyx != 0) { for (l = 0; l < 4; ++l) { pcs[l].p = node_pts + offset[l]; pcs[l].comp[0] = comp[offset[l]]; pcs[l].s[0] = states[offset[l]]; pcs[l].comp[1] = NO_COMP; pcs[l].s[1] = NULL; } num_pcs_els_in_blk(blk_el0) = 4; add_intfc_blk_pcs2d(*nbyx,*byx,*cyx,blk_el0, hash_table,h_size); num_pcs += num_pcs_els_in_blk(blk_el0); } else { set_bilinear_blk_el0(blk_el0); pcs[0].p = node_pts; pcs[0].comp[0] = comp[0]; pcs[0].s[0] = states[0]; pcs[0].comp[1] = NO_COMP; pcs[0].s[1] = NULL; ++num_pcs; } ++node_pts; ++states; ++comp; ++blk_el0; } ++node_pts; ++states; ++comp; /* if (iy < ymax -1) { node_pts += xmax+1; states += xmax+1; comp += xmax+1; } */ } free(hash_table); #if defined(DEBUG_TRI_GRID) if (debugging("collect_pcs_in_mesh2d")) { int i; int icoords[MAXD]; float *coords; (void) printf("\t\t PRINTING PCS:\n"); for (i=0;i<ntg->n_pcs ;++i) { coords = (float*)(ntg->pcs[i].p); if (rect_in_which(coords,icoords,&(ntg->rect_grid))) { (void) printf("%d: icoords[%d,%d] ", i,icoords[0],icoords[1]); (void) printf("comp[0] = %d, comp[1] = %d \n", ntg->pcs[i].comp[0],ntg->pcs[i].comp[1]); } else { (void) printf("rect_in_which failed for this pcs\n"); (void) printf("%d: comp[0] = %d, comp[1] = %d \n", i,ntg->pcs[i].comp[0],ntg->pcs[i].comp[1]); } (void) printf("(%g,%g)\n", coords[0],coords[1]); } (void) printf("\t\n\n PRINTING BLK_EL0:\n"); for (icoords[0] = 0; icoords[0] < ntg->rect_grid.gmax[0]; ++icoords[0]) { for (icoords[1] = 0; icoords[1] < ntg->rect_grid.gmax[1]; ++icoords[1]) { blk_el0 = &Blk_el0(icoords,ntg); (void) printf("blk_el0[%d,%d] = %llu\n", icoords[0],icoords[1], ptr2ull(blk_el0)); (void) printf("num_pcs_els_in_blk = %d\n", num_pcs_els_in_blk(blk_el0)); for (i = 0; i < num_pcs_els_in_blk(blk_el0); ++i) { (void) printf("%d: (%g, %g)\n",i, Coords(blk_el0_pcs_els(blk_el0)[i].p)[0], Coords(blk_el0_pcs_els(blk_el0)[i].p)[1]); } } } } if (debugging("tri_grid")) { (void) printf("Grid interface AFTER collect_pcs_in_mesh2d()\n"); print_interface(ntg->grid_intfc); } #endif /* defined(DEBUG_TRI_GRID) */ return GOOD_STEP; } /*end collect_pcs_in_mesh2d*/
EXPORT void init_topological_grid( RECT_GRID *top_grid, const RECT_GRID *r_grid) { char *c; const char *dimname; static const char *blanks = " \t"; const char *fmt1,*fmt2,*fmt3; char *message; double cor_fac; int len; int i, dim = r_grid->dim; int n_ints, n_floats, gmax[3]; static const char *dimnames[] = {"one", "two", "three"}; char s[Gets_BUF_SIZE]; dimname = dimnames[dim-1]; copy_rect_grid(top_grid,r_grid); (void) adjust_top_grid_for_square(top_grid,r_grid); for (i = 0; i < dim; ++i) gmax[i] = top_grid->gmax[i]; fmt1 = "The topological grid is a grid used for the construction of " "the tracked front topology. " "It is constrained to be a square grid. " "You specify the grid in one of two ways. " "If you enter a single number, it will be used as a " "coarseness factor for the topological grid relative to the " "computational grid entered above. " "In this case the length of a topological grid block cell " "side is the nearest allowable multiple of the shortest side " "of the computational grid by the coarseness factor. "; fmt2 = "Otherwise the code will read the %s integers input for the " "number of grid cells in each coordinate direction of the " "topological grid. " "If your input values do not yield a square grid they will " "be corrected to produce a square grid. " "This correction will attempt to produce values close to those " "input, but if the input values are highly rectangular, " "the resulting values may differ considerably from those " "entered. "; fmt3 = "The default for this input option is the nearest " "square grid that matches the computational grid. " "Generally the topological grid is coarser than the " "computational grid. " "Larger coarseness factors yield coarser grids, a value one " "gives the nearest square grid to the computational grid.\n"; len = 500; uni_array(&message,len,CHAR); (void) sprintf(message,fmt1,dimname); screen_print_long_string(message); (void) sprintf(message,fmt2,dimname); screen_print_long_string(message); (void) sprintf(message,fmt3,dimname); screen_print_long_string(message); free(message); message = NULL; if (dim == 1) { screen_print_long_string( "\tBe sure to use a decimal point to " "indicate a floating point number " "if you choose to input a coarseness factor.\n"); } screen("Enter your choice (cor_fac, %s integer%s, or return)\n", dimname,(dim > 1) ? "s" : ""); screen("\t(defaults are"); for (i = 0; i < dim; ++i) screen(" %d",gmax[i]); screen("): "); (void) Gets(s); if (s[0] != '\0') { n_floats = sscan_float(s,&cor_fac); if (dim == 1) { /* * To distinguish whether the input is a * coarseness factor search for a "." in s */ if (strstr(s,".") == NULL) { n_floats = 0; cor_fac = ERROR_FLOAT; } } for (n_ints = 0, c = strtok(s,blanks); c != NULL; c = strtok(NULL,blanks), ++n_ints) { (void) sscanf(c,"%d",gmax+n_ints%dim); } if (n_ints == 2*dim) { /* There is a secret undocumemted option here. * Enter the topological grid sizes twice and the * code will use these values whether they give a square grid * or not. This makes nearest_interface_point() * invalid, but you can get what you ask for. */ set_rect_grid(top_grid->L,top_grid->U,top_grid->L,top_grid->U, NOBUF,NOBUF,gmax,dim,&r_grid->Remap,top_grid); } else if ((n_ints == 1) && (n_floats > 0)) { int imin; imin = 0; for (i = 1; i < dim; ++i) if (r_grid->h[i] < r_grid->h[imin]) imin = i; top_grid->gmax[imin] = irint(r_grid->gmax[imin]/cor_fac); if (top_grid->gmax[imin] <= 0) top_grid->gmax[imin] = 1; top_grid->h[imin] = (top_grid->U[imin] - top_grid->L[imin]) / top_grid->gmax[imin]; for (i = 0; i < dim; ++i) { double tmp_gmax; if (i == imin) continue; tmp_gmax = (top_grid->U[i] - top_grid->L[i]) / top_grid->h[imin]; top_grid->gmax[i] = (int)(tmp_gmax); } (void) adjust_top_grid_for_square(top_grid,r_grid); set_rect_grid(top_grid->L,top_grid->U,top_grid->L,top_grid->U, NOBUF,NOBUF,top_grid->gmax,dim,&r_grid->Remap,top_grid); } else if (n_ints == dim) { for (i = 0; i < dim; ++i) top_grid->gmax[i] = gmax[i]; (void) adjust_top_grid_for_square(top_grid,r_grid); set_rect_grid(top_grid->L,top_grid->U,top_grid->L,top_grid->U, NOBUF,NOBUF,top_grid->gmax,dim,&r_grid->Remap,top_grid); } else { screen("ERROR in init_topological_grid(), " "invalid input of topogical grid mesh\n"); clean_up(ERROR); } } screen("The topological mesh used is "); for (i = 0; i < dim; ++i) screen(" %d",top_grid->gmax[i]); screen("\n\n"); } /*end init_topological_grid*/