void read_case(char *filename, int *n_variables, int *n_nodes, struct NODE **node, int *n_faces, struct FACE **face, int *n_cells, struct CELL **cell, int *n_zones, struct ZONE **zone) { int i; //open the file FILE *file = fopen(filename,"r"); exit_if_false(file != NULL, "opening the case file"); //number of variables exit_if_false(fread(n_variables, sizeof(int), 1, file) == 1, "reading the number of variables"); //numbers of elements exit_if_false(fread(n_nodes, sizeof(int), 1, file) == 1, "reading the number of nodes"); exit_if_false((*node = nodes_new(*n_nodes,NULL)) != NULL,"allocating nodes"); exit_if_false(fread(n_faces, sizeof(int), 1, file) == 1, "reading the number of faces"); exit_if_false((*face = faces_new(*n_faces,NULL)) != NULL,"allocating faces"); exit_if_false(fread(n_cells, sizeof(int), 1, file) == 1, "reading the number of cells"); exit_if_false((*cell = cells_new(*n_cells,NULL)) != NULL,"allocating cells"); exit_if_false(fread(n_zones, sizeof(int), 1, file) == 1, "reading the number of zones"); exit_if_false((*zone = zones_new(*n_zones,NULL)) != NULL,"allocating zones"); //elements for(i = 0; i < *n_nodes; i ++) node_case_get(file, &(*node)[i]); for(i = 0; i < *n_faces; i ++) face_case_get(file, *node, &(*face)[i], *cell, *zone); for(i = 0; i < *n_cells; i ++) cell_case_get(file, *n_variables, *face, &(*cell)[i], *zone); for(i = 0; i < *n_zones; i ++) zone_case_get(file, &(*zone)[i]); //clean up fclose(file); }
void zone_case_get(FILE *file, struct ZONE *zone) { exit_if_false(fread(&(zone->location), sizeof(char), 1, file) == 1,"reading the zone location"); exit_if_false(fread(&(zone->variable), sizeof(int), 1, file) == 1,"reading the zone variable"); exit_if_false(fread(zone->condition, sizeof(char), MAX_CONDITION_LENGTH, file) == MAX_CONDITION_LENGTH,"reading the zone condition"); exit_if_false(fread(&(zone->value), sizeof(double), 1, file) == 1,"reading the zone value"); }
void read_data(char *filename, double *time, int n_data, double *data) { FILE *file = fopen(filename,"r"); exit_if_false(file != NULL,"opening data file"); exit_if_false(fread(time, sizeof(double), 1, file) == 1,"reading the time"); exit_if_false(fread(data, sizeof(double), n_data, file) == n_data,"reading the data"); fclose(file); }
void generate_timed_named_filename(char *filename, char *basename, double time, char *name) { char *sub[2]; sub[0] = strchr(basename, '?'); exit_if_false(sub[0] != NULL,"first substitute character \"?\" not found in data basename"); sub[1] = strchr(sub[0] + 1, '?'); exit_if_false(sub[1] != NULL,"second substitute character \"?\" not found in data basename"); union { int integer; float real; } number; number.real = time; *sub[0] = *sub[1] = '\0'; sprintf(filename, "%s%i%s%s%s", basename, number.integer, sub[0] + 1, name, sub[1] + 1); *sub[0] = *sub[1] = '?'; }
void write_data(char *basename, double time, int n_data, double *data) { char *filename; exit_if_false(allocate_character_vector(&filename, MAX_STRING_LENGTH),"allocating filename"); generate_timed_filename(filename, basename, time); FILE *file = fopen(filename,"w"); exit_if_false(file != NULL,"opening data file"); exit_if_false(fwrite(&time, sizeof(double), 1, file) == 1,"writing the time"); exit_if_false(fwrite(data, sizeof(double), n_data, file) == n_data,"writing the data"); free_vector(filename); fclose(file); }
int update_element_integration(int n_variables_old, int n_variables, int *variable_order_old, int *variable_order, int n_elements, struct ELEMENT *element) { int e, f, h, i, j, k, o, v; int max_variable_order = 0, max_variable_order_old = 0; for(v = 0; v < n_variables; v ++) max_variable_order = MAX(max_variable_order,variable_order[v]); for(v = 0; v < n_variables_old; v ++) max_variable_order_old = MAX(max_variable_order_old,variable_order_old[v]); if(max_variable_order == max_variable_order_old) return 0; int n_hammer = ORDER_TO_N_HAMMER(max_variable_order); double dx[2][2], size; for(e = 0; e < n_elements; e ++) { exit_if_false(element[e].X = allocate_element_x(&element[e],n_hammer*(element[e].n_faces - 2)),"allocating element X"); exit_if_false(element[e].W = allocate_element_w(&element[e],n_hammer*(element[e].n_faces - 2)),"allocating element W"); f = 0; for(i = 0; i < element[e].n_faces - 2; i ++) { // triangulate f ++; while(element[e].face[f]->node[0] == element[e].face[0]->node[0] || element[e].face[f]->node[1] == element[e].face[0]->node[0]) f ++; o = element[e].face[f]->border[0] == &element[e]; for(j = 0; j < 2; j ++) for(k = 0; k < 2; k ++) dx[j][k] = element[e].face[f]->node[j == o]->x[k] - element[e].face[0]->node[0]->x[k]; // hammer locations and weights size = dx[0][0]*dx[1][1] - dx[0][1]*dx[1][0]; for(h = 0; h < n_hammer; h ++) { for(j = 0; j < 2; j ++) { element[e].X[j][i*n_hammer+h] = element[e].face[0]->node[0]->x[j] + hammer_locations[n_hammer-1][0][h]*dx[0][j] + hammer_locations[n_hammer-1][1][h]*dx[1][j]; } element[e].W[i*n_hammer+h] = hammer_weights[n_hammer-1][h] * size; } } } return n_elements; }
void write_gnuplot(char *basename, double time, int n_variables, char **variable_name, int n_unknowns, int *unknown_to_id, double *x, int n_faces, struct FACE *face, int n_cells, struct CELL *cell, int n_zones, struct ZONE *zone) { int i, j, u, v, z; int n_polygon; double ***polygon; exit_if_false(allocate_double_pointer_matrix(&polygon,MAX(MAX_CELL_FACES,4),2),"allocating polygon memory"); char *filename; exit_if_false(allocate_character_vector(&filename, MAX_STRING_LENGTH),"allocating filename"); FILE **file = (FILE **)malloc(n_variables * sizeof(FILE *)); exit_if_false(file != NULL,"allocating files"); for(v = 0; v < n_variables; v ++) { generate_timed_named_filename(filename, basename, time, variable_name[v]); file[v] = fopen(filename,"w"); exit_if_false(file[v] != NULL,"opening file"); } for(u = 0; u < n_unknowns; u ++) { i = ID_TO_INDEX(unknown_to_id[u]); z = ID_TO_ZONE(unknown_to_id[u]); n_polygon = generate_control_volume_polygon(polygon, i, zone[z].location, face, cell); v = zone[z].variable; for(j = 1; j < n_polygon; j ++) { fprintf(file[v],"%+.10e %+.10e %+.10e\n",polygon[0][0][0],polygon[0][0][1],x[u]); fprintf(file[v],"%+.10e %+.10e %+.10e\n\n",polygon[j][0][0],polygon[j][0][1],x[u]); fprintf(file[v],"%+.10e %+.10e %+.10e\n",polygon[0][0][0],polygon[0][0][1],x[u]); fprintf(file[v],"%+.10e %+.10e %+.10e\n\n\n",polygon[j][1][0],polygon[j][1][1],x[u]); } } for(v = 0; v < n_variables; v ++) fclose(file[v]); free_matrix((void **)polygon); free_vector(filename); free(file); }
void cell_geometry_get(FILE *file, struct CELL *cell, struct FACE *face) { int i; //temporary storage int *index, count, offset; char *line, *temp; index = (int *)malloc(MAX_CELL_FACES * sizeof(int)); line = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); temp = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); exit_if_false(index != NULL && line != NULL && temp != NULL ,"allocating temporary storage"); //read the line exit_if_false(fgets(line, MAX_STRING_LENGTH, file) != NULL, "reading a cell line"); //eat up whitespace and newlines for(i = strlen(line)-1; i >= 0; i --) if(line[i] != ' ' && line[i] != '\n') break; line[i+1] = '\0'; //sequentially read the integers on the line count = offset = 0; while(offset < strlen(line)) { sscanf(&line[offset],"%s",temp); sscanf(temp,"%i",&index[count]); count ++; offset += strlen(temp) + 1; while(line[offset] == ' ') offset ++; } //number of faces cell->n_faces = count; //allocate the faces exit_if_false(cell_face_new(cell),"allocating cell faces"); //face pointers for(i = 0; i < count; i ++) cell->face[i] = &face[index[i]]; //clean up free(index); free(line); free(temp); }
void face_geometry_get(FILE *file, struct FACE *face, struct NODE *node) { int i; //temporary storage int *index, count, offset; char *line, *temp; index = (int *)malloc(MAX_FACE_NODES * sizeof(int)); line = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); temp = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); exit_if_false(index != NULL && line != NULL && temp != NULL ,"allocating temporary storage"); //read a line exit_if_false(fgets(line, MAX_STRING_LENGTH, file) != NULL, "reading a face line"); //strip newlines and whitespace off the end of the line for(i = strlen(line)-1; i >= 0; i --) if(line[i] != ' ' && line[i] != '\n') break; line[i+1] = '\0'; //sequentially read the integers on the line count = offset = 0; while(offset < strlen(line)) { sscanf(&line[offset],"%s",temp); sscanf(temp,"%i",&index[count]); count ++; offset += strlen(temp) + 1; while(line[offset] == ' ') offset ++; } //number of faces face->n_nodes = count; //allocate the faces exit_if_false(face_node_new(face),"allocating face nodes"); //node pointers for(i = 0; i < count; i ++) face->node[i] = &node[index[i]]; //clean up free(index); free(line); free(temp); }
void generate_timed_filename(char *filename, char *basename, double time) { char *sub = strchr(basename, '?'); exit_if_false(sub != NULL,"substitute character \"?\" not found in data basename"); union { int integer; float real; } number; number.real = time; *sub = '\0'; sprintf(filename, "%s%i%s", basename, number.integer, sub + 1); *sub = '?'; }
void write_case(char *filename, int n_variables, int n_nodes, struct NODE *node, int n_faces, struct FACE *face, int n_cells, struct CELL *cell, int n_zones, struct ZONE *zone) { int i; //open the file FILE *file = fopen(filename,"w"); exit_if_false(file != NULL, "opening the case file"); //number of variables exit_if_false(fwrite(&n_variables, sizeof(int), 1, file) == 1, "writing the number of variables"); //numbers of elements exit_if_false(fwrite(&n_nodes, sizeof(int), 1, file) == 1, "writing the number of nodes"); exit_if_false(fwrite(&n_faces, sizeof(int), 1, file) == 1, "writing the number of faces"); exit_if_false(fwrite(&n_cells, sizeof(int), 1, file) == 1, "writing the number of cells"); exit_if_false(fwrite(&n_zones, sizeof(int), 1, file) == 1, "writing the number of zones"); //elements for(i = 0; i < n_nodes; i ++) node_case_write(file, &node[i]); for(i = 0; i < n_faces; i ++) face_case_write(file, node, &face[i], cell, zone); for(i = 0; i < n_cells; i ++) cell_case_write(file, n_variables, face, &cell[i], zone); for(i = 0; i < n_zones; i ++) zone_case_write(file, &zone[i]); //clean up fclose(file); }
void face_case_write(FILE *file, struct NODE *node, struct FACE *face, struct CELL *cell, struct ZONE *zone) { int i, *index; index = (int *)malloc(MAX_INDICES * sizeof(int)); exit_if_false(index != NULL,"allocating temporary storage"); exit_if_false(fwrite(&(face->n_nodes), sizeof(int), 1, file) == 1, "writing the number of face nodes"); for(i = 0; i < face->n_nodes; i ++) index[i] = (int)(face->node[i] - &node[0]); exit_if_false(fwrite(index, sizeof(int), face->n_nodes, file) == face->n_nodes, "writing the face nodes"); exit_if_false(fwrite(&(face->area), sizeof(double), 1, file) == 1, "writing the face area"); exit_if_false(fwrite(face->centroid, sizeof(double), 2, file) == 2, "writing the face centroid"); exit_if_false(fwrite(&(face->n_borders), sizeof(int), 1, file) == 1, "writing the number of face borders"); for(i = 0; i < face->n_borders; i ++) index[i] = (int)(face->border[i] - &cell[0]); exit_if_false(fwrite(index, sizeof(int), face->n_borders, file) == face->n_borders, "writing the face borders"); exit_if_false(fwrite(face->oriented, sizeof(int), face->n_borders, file) == face->n_borders, "writing the face orientations"); exit_if_false(fwrite(&(face->n_zones), sizeof(int), 1, file) == 1, "writing the number of face zones"); for(i = 0; i < face->n_zones; i ++) index[i] = (int)(face->zone[i] - &zone[0]); exit_if_false(fwrite(index, sizeof(int), face->n_zones, file) == face->n_zones, "writing the face zones"); free(index); }
int update_face_integration(int n_variables_old, int n_variables, int *variable_order_old, int *variable_order, int n_faces, struct FACE *face) { int f, g, i, v; int max_variable_order = 0, max_variable_order_old = 0; for(v = 0; v < n_variables; v ++) max_variable_order = MAX(max_variable_order,variable_order[v]); for(v = 0; v < n_variables_old; v ++) max_variable_order_old = MAX(max_variable_order_old,variable_order_old[v]); if(max_variable_order == max_variable_order_old) return 0; int n_gauss = ORDER_TO_N_GAUSS(max_variable_order); double dx[2], size; for(f = 0; f < n_faces; f ++) { exit_if_false(face[f].X = allocate_face_x(&face[f], ORDER_TO_N_GAUSS(max_variable_order)),"allocating face X"); exit_if_false(face[f].W = allocate_face_w(&face[f], ORDER_TO_N_GAUSS(max_variable_order)),"allocating face W"); for(i = 0; i < 2; i ++) dx[i] = face[f].node[1]->x[i] - face[f].node[0]->x[i]; size = sqrt(dx[0]*dx[0] + dx[1]*dx[1]); for(g = 0; g < n_gauss; g ++) { for(i = 0; i < 2; i ++) { face[f].X[i][g] = 0.5 * (1.0 - gauss_locations[n_gauss-1][g]) * face[f].node[0]->x[i] + 0.5 * (1.0 + gauss_locations[n_gauss-1][g]) * face[f].node[1]->x[i]; } face[f].W[g] = 0.5 * size * gauss_weights[n_gauss-1][g]; } } return n_faces; }
void write_vtk(char *basename, double time, int n_variables, char **variable_name, int n_unknowns, int *unknown_to_id, double *x, int n_nodes, struct NODE *node, int n_faces, struct FACE *face, int n_cells, struct CELL *cell, int n_zones, struct ZONE *zone) { char *filename; exit_if_false(allocate_character_vector(&filename, MAX_STRING_LENGTH),"allocating filename"); FILE *file; int *point_used, *point_index; exit_if_false(allocate_integer_vector(&point_used,n_nodes + n_cells),"allocating point usage array"); exit_if_false(allocate_integer_vector(&point_index,n_nodes + n_cells),"allocating point index array"); int n_points, n_elements; int i, j, u, v, z, offset; for(v = 0; v < n_variables; v ++) { generate_timed_named_filename(filename, basename, time, variable_name[v]); file = fopen(filename,"w"); exit_if_false(file != NULL,"opening file"); fprintf(file,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n<UnstructuredGrid>\n"); for(i = 0; i < n_nodes + n_cells; i ++) point_used[i] = 0; n_elements = 0; for(u = 0; u < n_unknowns; u ++) { i = ID_TO_INDEX(unknown_to_id[u]); z = ID_TO_ZONE(unknown_to_id[u]); if(zone[z].variable == v) { n_elements ++; if(zone[z].location == 'f') { point_used[(int)(face[i].node[0] - &node[0])] = 1; point_used[(int)(face[i].node[face[i].n_nodes - 1] - &node[0])] = 1; for(j = 0; j < face[i].n_borders; j ++) point_used[n_nodes + (int)(face[i].border[j] - &cell[0])] = 1; } else if(zone[z].location == 'c') { for(j = 0; j < cell[i].n_faces; j ++) { point_used[(int)(cell[i].face[j]->node[0] - &node[0])] = 1; point_used[(int)(face[i].node[cell[i].face[j]->n_nodes - 1] - &node[0])] = 1; } } else exit_if_false(0,"recognising the location"); } } n_points = 0; for(i = 0; i < n_nodes + n_cells; i ++) if(point_used[i]) point_index[n_points ++] = i; fprintf(file,"<Piece NumberOfPoints=\"%i\" NumberOfCells=\"%i\">\n", n_points, n_elements); fprintf(file,"<CellData>\n"); fprintf(file,"<DataArray Name=\"%s\" type=\"Float64\" format=\"ascii\">\n",variable_name[v]); for(u = 0; u < n_unknowns; u ++) if(zone[ID_TO_ZONE(unknown_to_id[u])].variable == v) fprintf(file," %+.10e",x[u]); fprintf(file,"\n</DataArray>\n"); fprintf(file,"</CellData>\n"); fprintf(file,"<Points>\n"); fprintf(file,"<DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">\n"); for(i = 0; i < n_nodes; i ++) if(point_used[i]) fprintf(file," %+.10e %+.10e %+.10e",node[i].x[0],node[i].x[1],0.0); for(i = 0; i < n_cells; i ++) if(point_used[n_nodes + i]) fprintf(file," %+.10e %+.10e %+.10e",cell[i].centroid[0],cell[i].centroid[1],0.0); fprintf(file,"\n</DataArray>\n"); fprintf(file,"</Points>\n"); fprintf(file,"<Cells>\n"); fprintf(file,"<DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">\n"); for(u = 0; u < n_unknowns; u ++) { i = ID_TO_INDEX(unknown_to_id[u]); z = ID_TO_ZONE(unknown_to_id[u]); if(zone[z].variable == v) { if(zone[z].location == 'f') { fprintf(file," %i",point_index[(int)((face[i].oriented[0] ? face[i].node[1] : face[i].node[0]) - &node[0])]); fprintf(file," %i",point_index[n_nodes + (int)(face[i].border[0] - &cell[0])]); fprintf(file," %i",point_index[(int)((face[i].oriented[0] ? face[i].node[0] : face[i].node[1]) - &node[0])]); if(face[i].n_borders == 2) fprintf(file," %i",point_index[n_nodes + (int)(face[i].border[1] - &cell[0])]); } else if(zone[z].location == 'c') { for(j = 0; j < cell[i].n_faces; j ++) { fprintf(file," %i",point_index[(int)(cell[i].face[j]->node[!cell[i].oriented[j]] - &node[0])]); } } else exit_if_false(0,"recognising the location"); } } fprintf(file,"\n</DataArray>\n"); fprintf(file,"<DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">\n"); offset = 0; for(u = 0; u < n_unknowns; u ++) { i = ID_TO_INDEX(unknown_to_id[u]); z = ID_TO_ZONE(unknown_to_id[u]); if(zone[z].variable == v) { if(zone[z].location == 'f') { offset += 2 + face[i].n_borders; } else if(zone[z].location == 'c') { offset += cell[i].n_faces; } else exit_if_false(0,"recognising the location"); fprintf(file," %i",offset); } } fprintf(file,"\n</DataArray>\n"); fprintf(file,"<DataArray type=\"Int32\" Name=\"types\" format=\"ascii\">\n"); for(i = 0; i < n_elements; i ++) fprintf(file," %i",7); fprintf(file,"\n</DataArray>\n"); fprintf(file,"</Cells>"); fprintf(file,"\n</Piece>\n</UnstructuredGrid>\n</VTKFile>"); fclose(file); } free_vector(filename); free_vector(point_used); free_vector(point_index); }
int update_element_numerics(int n_variables_old, int n_variables, int *variable_order_old, int *variable_order, int n_elements, struct ELEMENT *element) { int e, i, j, v; // old and new maximum variable orders int max_variable_order = 0, max_variable_order_old = 0; for(v = 0; v < n_variables; v ++) max_variable_order = MAX(max_variable_order,variable_order[v]); for(v = 0; v < n_variables_old; v ++) max_variable_order_old = MAX(max_variable_order_old,variable_order_old[v]); // what needs updating int max_update = max_variable_order != max_variable_order_old, any_update = n_variables_old < n_variables; for(v = 0; v < MIN(n_variables_old,n_variables); v ++) any_update = any_update || (variable_order_old[v] != variable_order[v]); if(!any_update) return 0; // numbers of basis functions int *n_basis, max_n_basis = ORDER_TO_N_BASIS(max_variable_order); exit_if_false(n_basis = (int *)malloc(n_variables * sizeof(int)),"allocating n_basis"); for(v = 0; v < n_variables; v ++) n_basis[v] = ORDER_TO_N_BASIS(variable_order[v]); // numbers of points int n_gauss = ORDER_TO_N_GAUSS(max_variable_order), n_hammer = ORDER_TO_N_HAMMER(max_variable_order), n_points; // no differential int no_differential[2] = {0,0}; // working matrices int lds = max_n_basis, ldm = max_n_basis, ldd = max_n_basis, lda = max_n_basis; int sizem = max_n_basis*max_n_basis; double **S = allocate_double_matrix(NULL,(MAX_ELEMENT_N_FACES-2)*n_hammer,lds); double **M = allocate_double_matrix(NULL,max_n_basis,ldm); double **D = allocate_double_matrix(NULL,max_n_basis,ldd); double **A = allocate_double_matrix(NULL,max_n_basis,lda); double **X = allocate_double_matrix(NULL,2,MAX_ELEMENT_N_FACES); exit_if_false(S != NULL && M != NULL && A != NULL && X != NULL,"allocating working matrices"); // lapack and blas int info, *pivot = (int *)malloc((max_n_basis + 2) * sizeof(int)); exit_if_false(pivot != NULL,"allocating pivot"); char trans[2] = "NT"; int int_1 = 1; double dbl_0 = 0.0, dbl_1 = 1.0; for(e = 0; e < n_elements; e ++) { n_points = n_hammer*(element[e].n_faces - 2); if(max_update) { // interior matrices exit_if_false(element[e].P = allocate_element_p(&element[e],max_n_basis,n_points),"allocating element P"); for(i = 0; i < max_n_basis; i ++) for(j = 0; j < max_n_basis; j ++) basis(n_points,element[e].P[i][j],element[e].X,element[e].centre,element[e].size,j,taylor_powers[i]); // face matrices exit_if_false(element[e].Q = allocate_element_q(&element[e],max_n_basis,n_gauss),"allocating element Q"); for(i = 0; i < element[e].n_faces; i ++) for(j = 0; j < max_n_basis; j ++) basis(n_gauss,element[e].Q[i][j],element[e].face[i]->X,element[e].centre,element[e].size,j,no_differential); // corner matrix exit_if_false(element[e].V = allocate_element_v(&element[e],max_n_basis),"allocating element V"); for(i = 0; i < element[e].n_faces; i ++) for(j = 0; j < 2; j ++) X[j][i] = element[e].face[i]->node[element[e].face[i]->border[0] != &element[e]]->x[j]; for(i = 0; i < max_n_basis; i ++) basis(element[e].n_faces,element[e].V[i],X,element[e].centre,element[e].size,i,no_differential); } // mass matrix for(i = 0; i < max_n_basis; i ++) dcopy_(&n_points,element[e].P[powers_taylor[0][0]][i],&int_1,&S[0][i],&lds); for(i = 0; i < n_points; i ++) dscal_(&max_n_basis,&element[e].W[i],S[i],&int_1); dgemm_(&trans[0],&trans[0],&max_n_basis,&max_n_basis,&n_points,&dbl_1,S[0],&lds,element[e].P[powers_taylor[0][0]][0],&n_points,&dbl_0,M[0],&ldm); // initialise matrices exit_if_false(element[e].I = allocate_element_i(&element[e],n_variables,n_basis,n_points),"allocating element I"); for(v = 0; v < n_variables; v ++) { if(!max_update && n_variables_old > v) if(variable_order_old[v] == variable_order[v]) continue; for(i = 0; i < n_basis[v]; i ++) dcopy_(&n_points,&S[0][i],&lds,&element[e].I[v][0][i],&n_basis[v]); dcopy_(&sizem,M[0],&int_1,A[0],&int_1); dgesv_(&n_basis[v],&n_points,A[0],&lda,pivot,element[e].I[v][0],&n_basis[v],&info); } // limiting matrices if(max_n_basis > 1) { // diffusion matrix for(i = 0; i < max_n_basis; i ++) dcopy_(&n_points,element[e].P[powers_taylor[1][0]][i],&int_1,&S[0][i],&lds); for(i = 0; i < n_points; i ++) dscal_(&max_n_basis,&element[e].W[i],S[i],&int_1); dgemm_(&trans[0],&trans[0],&max_n_basis,&max_n_basis,&n_points,&dbl_1,S[0],&lds,element[e].P[powers_taylor[1][0]][0],&n_points,&dbl_0,D[0],&ldd); for(i = 0; i < max_n_basis; i ++) dcopy_(&n_points,element[e].P[powers_taylor[0][1]][i],&int_1,&S[0][i],&lds); for(i = 0; i < n_points; i ++) dscal_(&max_n_basis,&element[e].W[i],S[i],&int_1); dgemm_(&trans[0],&trans[0],&max_n_basis,&max_n_basis,&n_points,&dbl_1,S[0],&lds,element[e].P[powers_taylor[0][1]][0],&n_points,&dbl_1,D[0],&ldd); } // limiting matrices exit_if_false(element[e].L = allocate_element_l(&element[e],n_variables,n_basis),"allocating element L"); for(v = 0; v < n_variables; v ++) { if(n_variables_old > v) if(variable_order_old[v] == variable_order[v]) continue; if(variable_order[v] == 1) continue; dcopy_(&sizem,M[0],&int_1,A[0],&int_1); for(i = 0; i < n_basis[v]; i ++) dcopy_(&n_basis[v],D[i],&int_1,element[e].L[v][i],&int_1); dgesv_(&n_basis[v],&n_basis[v],A[0],&lda,pivot,element[e].L[v][0],&n_basis[v],&info); } } free(n_basis); destroy_matrix((void *)S); destroy_matrix((void *)M); destroy_matrix((void *)D); destroy_matrix((void *)A); destroy_matrix((void *)X); free(pivot); return n_elements; }
int update_face_numerics(int n_variables_old, int n_variables, int *variable_order_old, int *variable_order, int n_faces, struct FACE *face, int n_boundaries_old, struct BOUNDARY *boundary_old) { int a, b, f, g, h, i, j, v; // lists of old face boundaries int **face_n_boundaries_old = allocate_integer_matrix(NULL,n_faces,n_variables); exit_if_false(face_n_boundaries_old != NULL,"allocating face_n_boundaries_old"); for(i = 0; i < n_faces; i ++) for(j = 0; j < n_variables; j ++) face_n_boundaries_old[i][j] = 0; int ***face_boundary_old = allocate_integer_tensor(NULL,n_faces,n_variables,MAX_FACE_N_BOUNDARIES); exit_if_false(face_boundary_old != NULL,"allocating face_boundary_old"); for(b = 0; b < n_boundaries_old; b ++) { for(i = 0; i < boundary_old[b].n_faces; i ++) { v = boundary_old[b].variable; f = boundary_old[b].face[i] - &face[0]; face_boundary_old[f][v][face_n_boundaries_old[f][v]++] = b; } } // order and numbers of integration points and bases int max_variable_order = 0, max_variable_order_old = 0; for(v = 0; v < n_variables; v ++) max_variable_order = MAX(max_variable_order,variable_order[v]); for(v = 0; v < n_variables_old; v ++) max_variable_order_old = MAX(max_variable_order_old,variable_order_old[v]); int n_gauss = ORDER_TO_N_GAUSS(max_variable_order), n_hammer = ORDER_TO_N_HAMMER(max_variable_order); int n_points[MAX_FACE_N_BORDERS], sum_n_points[MAX_FACE_N_BORDERS+1]; int *n_basis = (int *)malloc(n_variables * sizeof(int)), max_n_basis = ORDER_TO_N_BASIS(max_variable_order); exit_if_false(n_basis != NULL,"allocating n_basis"); for(v = 0; v < n_variables; v ++) n_basis[v] = ORDER_TO_N_BASIS(variable_order[v]); // truth values for updating the interpolation on a face int *update = (int *)malloc(n_variables * sizeof(int)), max_update = max_variable_order != max_variable_order_old, any_update; exit_if_false(update != NULL,"allocating update"); // transformation double **R = allocate_double_matrix(NULL,2,2), det_R, **inv_R = allocate_double_matrix(NULL,2,2); exit_if_false(R != NULL,"allocating R"); exit_if_false(inv_R != NULL,"allocating inv_R"); double **T = allocate_double_matrix(NULL,max_n_basis,max_n_basis); exit_if_false(T != NULL,"allocating T"); // centre double **centre = allocate_double_matrix(NULL,2,1); exit_if_false(centre != NULL,"allocating centre"); // differentials int condition[2], none[2] = {0,0}; // integration locations in cartesian (x) and tranformed (y) coordinates double **x, **y; double **x_adj[MAX_FACE_N_BORDERS], ***y_adj; exit_if_false(y = allocate_double_matrix(NULL,2,n_gauss),"allocating y"); exit_if_false(y_adj = allocate_double_tensor(NULL,MAX_FACE_N_BORDERS,2,n_hammer*(MAX_ELEMENT_N_FACES-2)),"allocating y_adj"); // adjacent elements and boundaries to the face int n_adj, n_bnd; struct ELEMENT **adj; struct BOUNDARY **bnd; // interpolation problem sizes int n_adj_bases = MAX_FACE_N_BORDERS*max_n_basis; // number of adjacent bases int n_int_terms = MAX_FACE_N_BORDERS*max_n_basis + MAX_FACE_N_BOUNDARIES; // number of terms from which to interpolate int n_int_bases = MAX_FACE_N_BORDERS*max_n_basis + MAX_FACE_N_BOUNDARIES*max_variable_order; // number of interpolant bases // face basis taylor indices int *face_taylor = (int *)malloc(n_int_bases * sizeof(int)); exit_if_false(face_taylor != NULL,"allocating face_taylor"); // temporary matrices int ldp, lds, ldq; ldp = lds = ldq = MAX_FACE_N_BORDERS*(MAX_ELEMENT_N_FACES-2)*n_hammer; int sizep = n_adj_bases*ldp; double **P = allocate_double_matrix(NULL,n_adj_bases,ldp); double **S = allocate_double_matrix(NULL,n_adj_bases,lds); double **Q = allocate_double_matrix(NULL,n_int_bases,ldp); int lda, ldb; lda = ldb = n_int_bases; double **A = allocate_double_matrix(NULL,n_int_bases,n_int_bases); double **B = allocate_double_matrix(NULL,n_int_bases,n_int_bases); int ldf, ldd; ldf = ldd = n_gauss; double **F = allocate_double_matrix(NULL,n_int_bases,n_gauss); double ***D = allocate_double_tensor(NULL,max_n_basis,n_int_terms,n_gauss); int incd = n_int_terms*n_gauss, incq; exit_if_false(P != NULL && S != NULL && Q != NULL && A != NULL && B != NULL && F != NULL && D != NULL,"allocating working matrices"); // blas char trans[2] = "NT"; int i_one = 1; double d_one = 1.0, d_zero = 0.0; // lapack int info, *pivot = (int *)malloc(n_int_bases * sizeof(int)); exit_if_false(pivot != NULL,"allocating pivot"); int updated = 0; for(f = 0; f < n_faces; f ++) { // see if face interpolation needs updating any_update = 0; for(v = 0; v < n_variables; v ++) { update[v] = 0; if(v >= n_variables_old) update[v] = 1; else if(variable_order[v] != variable_order_old[v]) update[v] = 1; else if(face[f].n_boundaries[v] != face_n_boundaries_old[f][v]) update[v] = 1; else for(i = 0; i < face[f].n_boundaries[v]; i ++) for(j = 0; j < 2; j ++) if(face[f].boundary[v][i]->condition[j] != boundary_old[face_boundary_old[f][v][i]].condition[j]) update[v] = 1; any_update += update[v]; } if(!any_update) continue; // allocate matrices exit_if_false(face[f].Q = allocate_face_q(&face[f],n_variables,n_basis,n_gauss),"allocating face Q"); // rotation to face coordinates R[0][0] = + face[f].normal[0]; R[0][1] = + face[f].normal[1]; R[1][0] = - face[f].normal[1]; R[1][1] = + face[f].normal[0]; det_R = R[0][0]*R[1][1] - R[0][1]*R[1][0]; inv_R[0][0] = + R[1][1]/det_R; inv_R[0][1] = - R[0][1]/det_R; inv_R[1][0] = - R[1][0]/det_R; inv_R[1][1] = + R[0][0]/det_R; transformation_matrix(max_variable_order,T,inv_R); // face integration locations x = face[f].X; for(g = 0; g < n_gauss; g ++) { for(i = 0; i < 2; i ++) { y[i][g] = face[f].centre[i]; for(j = 0; j < 2; j ++) y[i][g] += R[i][j]*(x[j][g] - face[f].centre[j]); } } // numbers of element integration locations for(a = 0; a < face[f].n_borders; a ++) n_points[a] = n_hammer*(face[f].border[a]->n_faces-2); sum_n_points[0] = 0; for(a = 0; a < face[f].n_borders; a ++) sum_n_points[a+1] = sum_n_points[a] + n_points[a]; // adjacent element integration locations for(a = 0; a < face[f].n_borders; a ++) { x_adj[a] = face[f].border[a]->X; for(h = 0; h < n_points[a]; h ++) { for(i = 0; i < 2; i ++) { y_adj[a][i][h] = face[f].centre[i]; for(j = 0; j < 2; j ++) y_adj[a][i][h] += R[i][j]*(x_adj[a][j][h] - face[f].centre[j]); } } } // for all variables for(v = 0; v < n_variables; v ++) { if(!max_update && !update[v]) continue; n_adj = face[f].n_borders; adj = face[f].border; n_bnd = face[f].n_boundaries[v]; bnd = face[f].boundary[v]; n_adj_bases = n_adj*n_basis[v]; n_int_terms = n_adj_bases + n_bnd; // face basis indices n_int_bases = 0; for(i = 0; i < n_adj*variable_order[v] + n_bnd; i ++) for(j = 0; j < n_adj*variable_order[v] + n_bnd; j ++) if(i + n_adj*j < n_adj*variable_order[v] + n_bnd && j < variable_order[v]) face_taylor[n_int_bases ++] = powers_taylor[i][j]; exit_if_false(n_int_bases == n_adj_bases + n_bnd*variable_order[v],"mismatched number of interpolation unknowns"); // element bases at the integration locations for(i = 0; i < n_adj*n_basis[v]; i ++) for(j = 0; j < sum_n_points[n_adj]; j ++) P[i][j] = 0.0; for(a = 0; a < n_adj; a ++) for(i = 0; i < n_basis[v]; i ++) basis(n_points[a],&P[i+a*n_basis[v]][sum_n_points[a]],x_adj[a],adj[a]->centre,adj[a]->size,i,none); // face bases at the integration locations for(a = 0; a < n_adj; a ++) for(i = 0; i < n_int_bases; i ++) basis(n_points[a],&Q[i][sum_n_points[a]],y_adj[a],face[f].centre,face[f].size,face_taylor[i],none); // centre of face in form which can be passed to basis for(i = 0; i < 2; i ++) centre[i][0] = face[f].centre[i]; // integration matrix dcopy_(&sizep,P[0],&i_one,S[0],&i_one); for(a = 0; a < n_adj; a ++) for(i = 0; i < n_points[a]; i ++) dscal_(&n_basis[v],&adj[a]->W[i],&S[a*n_basis[v]][i+sum_n_points[a]],&lds); // weak interpolation system dgemm_(&trans[1],&trans[0],&n_adj_bases,&n_int_bases,&sum_n_points[n_adj],&d_one,S[0],&lds,Q[0],&ldq,&d_zero,A[0],&lda); // weak interpolation rhs dgemm_(&trans[1],&trans[0],&n_adj_bases,&n_adj_bases,&sum_n_points[n_adj],&d_one,S[0],&lds,P[0],&ldp,&d_zero,B[0],&ldb); // boundary conditions for(b = 0; b < n_bnd; b ++) { condition[0] = bnd[b]->condition[0]; for(i = 0; i < variable_order[v]; i ++) { condition[1] = bnd[b]->condition[1] + i; for(j = 0; j < n_int_bases; j ++) basis(1,&A[j][i+n_adj_bases],centre,face[f].centre,face[f].size,face_taylor[j],condition); } for(i = 0; i < variable_order[v]; i ++) for(j = 0; j < n_int_terms; j ++) B[j][i+n_adj_bases] = 0.0; for(i = 0; i < n_adj_bases; i ++) B[n_adj_bases+b][i] = 0.0; B[b+n_adj_bases][b*variable_order[v]+n_adj_bases] = 1.0; } // solve interpolation problem dgesv_(&n_int_bases,&n_int_terms,A[0],&lda,pivot,B[0],&ldb,&info); // interpolate values to the face integration locations for(i = 0; i < n_basis[v]; i ++) { for(j = 0; j < n_int_bases; j ++) basis(n_gauss,F[j],y,face[f].centre,face[f].size,face_taylor[j],taylor_powers[i]); dgemm_(&trans[0],&trans[0],&n_gauss,&n_int_terms,&n_int_bases,&d_one,F[0],&ldf,B[0],&ldb,&d_zero,D[i][0],&ldd); } // transform from face to cartesian coordinates incq = n_int_terms*n_gauss; for(i = 0; i < n_int_terms; i ++) for(j = 0; j < n_gauss; j ++) dgemv_(&trans[1],&n_basis[v],&n_basis[v],&d_one,T[0],&max_n_basis,&D[0][i][j],&incd,&d_zero,&face[f].Q[v][0][i][j],&incq); updated ++; } } // clean up destroy_matrix((void *)face_n_boundaries_old); destroy_tensor((void *)face_boundary_old); free(n_basis); free(update); destroy_matrix((void *)R); destroy_matrix((void *)inv_R); destroy_matrix((void *)T); destroy_matrix((void *)centre); destroy_matrix((void *)y); destroy_tensor((void *)y_adj); free(face_taylor); destroy_matrix((void *)P); destroy_matrix((void *)Q); destroy_matrix((void *)S); destroy_matrix((void *)A); destroy_matrix((void *)B); destroy_matrix((void *)F); destroy_tensor((void *)D); free(pivot); return updated; }
/** Runs the precision/recall test. * May log errors and even end the application in case of severe error. * @param params The program options. */ void eval_precision_recall( const program_options& params) { LOG(info) << "Loading class membership mappings..."; Vec1UInt membership_mappings; Vec1str cluster_files; exit_if_false( from_file( params.membership_mappings_file, membership_mappings), RETURN_CODE::IO_ERROR); exit_if_false( from_file( params.cluster_file_paths_file, cluster_files), RETURN_CODE::IO_ERROR); const uint n_features = static_cast<uint>(membership_mappings.size()); const uint n_clusters = static_cast<uint>(cluster_files.size()); LOG(info) << "# features: " << n_features; LOG(info) << "# clusters: " << n_clusters; LOG( info) << "calculating precision / recall ..."; // find best class for each cluster vector<cluster_info_t> class_mapping; for( uint i=0; i<n_clusters; ++i) { Vec1str current_cluster_image_paths; Vec1str current_cluster_real_image_classes; std::map<string, uint> class_votes; from_file( cluster_files[i], current_cluster_image_paths); for( auto it=current_cluster_image_paths.begin(); it!=current_cluster_image_paths.end(); ++it) { const string class_name = bfs::path(*it).parent_path().filename().string(); current_cluster_real_image_classes.push_back( class_name); const auto map_it = class_votes.find( class_name); if( map_it == class_votes.end()) class_votes[class_name] = 1; else map_it->second += 1; } const auto max_it = std::max_element( class_votes.begin(), class_votes.end(), []( const std::pair<string, int>& p, const std::pair<string, int>& q) { return p.second < q.second; }); if( max_it == class_votes.end()) { // ***cluster empty *** (yes, that can happen!) continue; } const string assigned_class( max_it->first); const uint n_retrieved_images( static_cast<uint>(current_cluster_image_paths.size())); uint false_positives(0); uint true_positives(0); uint false_negatives(0); // find true positives, false positives for( auto it=current_cluster_real_image_classes.begin(); it!=current_cluster_real_image_classes.end(); ++it) { const string& real_class = *it; if( assigned_class.compare( real_class) == 0) { ++true_positives; } else { ++false_positives; } } // false negatives std::stringstream folder_name; folder_name << params.image_db_directory << '/' << assigned_class; bfs::path folder_path( folder_name.str()); assert( bfs::exists( folder_path) && "the directory must exist."); uint n_relevant_images(0); for( bfs::directory_iterator it(folder_path); it!=bfs::directory_iterator(); ++it) { bfs::path p(*it); p.make_preferred(); if( !p.has_extension() || !is_image_filetype_supported( p.extension().string())) continue; ++n_relevant_images; Vec1str::iterator pos = std::find( current_cluster_image_paths.begin(), current_cluster_image_paths.end(), p.string()); if( pos == current_cluster_image_paths.end()) ++false_negatives; } assert( true_positives + false_negatives == n_relevant_images && "number of relevant images must be identical to the number of true positivies and false positives"); // *** found true positivies, false positives, false negatives for cluster i *** // calc precision/recall for each class const real precision = static_cast<real>(true_positives) / n_retrieved_images; const real recall = static_cast<real>(true_positives) / n_relevant_images; class_mapping.push_back( cluster_info_t( assigned_class, true_positives, false_positives, false_negatives, precision, recall)); } real avg_true_positives(0); real avg_false_positives(0); real avg_false_negatives(0); real avg_precision(0); real avg_recall(0); LOG( info) << "<class name> <true positives> <false positives> <false negatives> <precision> <recall>"; for( uint i=0; i<class_mapping.size(); ++i) { const cluster_info_t& ci = class_mapping[i]; const string& cluster_name = std::get<0>(ci); const uint true_positives = std::get<1>(ci); const uint false_positives = std::get<2>(ci); const uint false_negatives = std::get<3>(ci); const real precision = std::get<4>(ci); const real recall = std::get<5>(ci); avg_true_positives += true_positives; avg_false_positives += false_positives; avg_false_negatives += false_negatives; avg_precision += precision; avg_recall += recall; LOG(info) << cluster_name << " " << true_positives << " " << false_positives << " " << false_negatives << " " << precision << " " << recall; } avg_true_positives /= class_mapping.size(); avg_false_positives /= class_mapping.size(); avg_false_negatives /= class_mapping.size(); avg_precision /= class_mapping.size(); avg_recall /= class_mapping.size(); LOG(info) << "Average: <true positives> <false positives> <false negatives> <precision> <recall>"; LOG(info) << avg_true_positives << " " << avg_false_positives << " " << avg_false_negatives << " " << avg_precision << " " << avg_recall; LOG(info) << "Writing stats to file \"" << params.precision_recall_file << "\"..."; to_file( params.precision_recall_file, class_mapping); }
int main(int argc, char *argv[]) { exit_if_false(argc == 2,"exactly one input argument required"); timer_start(); //-------------------------------------------------------------------// // counters int i, n; // file paths char *input_file_path, *geometry_file_path, *case_file_path, *data_file_path, *data_numbered_file_path, *display_file_path, *display_numbered_file_path; exit_if_false(geometry_file_path = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating geometry file path"); exit_if_false(case_file_path = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating case file path"); exit_if_false(data_file_path = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating data file path"); exit_if_false(data_numbered_file_path = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating data numbered file path"); exit_if_false(display_file_path = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating display file path"); exit_if_false(display_numbered_file_path = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating display numbered file path"); // print string char *print; exit_if_false(print = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)),"allocating the print string"); // mesh structures int n_nodes, n_faces, n_elements, n_boundaries_old = 0, n_boundaries, n_terms; struct NODE *node; struct FACE *face; struct ELEMENT *element; struct BOUNDARY *boundary_old = NULL, *boundary; struct TERM *term; EXPRESSION *initial = NULL; // solution vectors int n_u_old = 0, n_u; double *u_old = NULL, *u; // files FILE *input_file, *case_file, *data_file, *geometry_file, *display_file; //-------------------------------------------------------------------// // opening the input file print_info("opening the input file %s",argv[1]); input_file_path = argv[1]; input_file = fopen(input_file_path,"r"); exit_if_false(input_file != NULL,"opening %s",input_file_path); // reading the case file path print_info("reading the case file path"); exit_if_false(fetch_value(input_file,"case_file_path",'s',case_file_path) == FETCH_SUCCESS,"reading case_file_path from %s",input_file_path); print_continue(case_file_path); // reading the number of variables, variable names and orders print_info("reading the variables"); int n_variables_old = 0, n_variables; exit_if_false(fetch_value(input_file,"number_of_variables",'i',&n_variables) == FETCH_SUCCESS,"reading number_of_variables from %s",input_file_path); int *variable_order_old = NULL, *variable_order = (int *)malloc(n_variables * sizeof(int)); exit_if_false(variable_order != NULL,"allocating orders"); exit_if_false(fetch_vector(input_file, "variable_order", 'i', n_variables, variable_order) == FETCH_SUCCESS,"reading variable_order from %s",input_file_path); char **variable_name = allocate_character_matrix(NULL,n_variables,MAX_STRING_LENGTH); exit_if_false(variable_name != NULL,"allocating variable names"); warn_if_false(fetch_vector(input_file,"variable_name",'s',n_variables,variable_name) == FETCH_SUCCESS,"reading variable_name from %s",input_file_path); for(i = 0; i < n_variables; i ++) print_continue("%s order %i",variable_name[i],variable_order[i]); // reading the number of inner and outer iterations to perform print_info("reading the numbers of iterations"); int outer_iteration = 0, inner_iteration; int n_outer_iterations, n_inner_iterations, data_n_outer_iterations, display_n_outer_iterations; exit_if_false(fetch_value(input_file,"number_of_inner_iterations",'i',&n_inner_iterations) == FETCH_SUCCESS,"reading number_of_inner_iterations from %s",input_file_path); exit_if_false(fetch_value(input_file,"number_of_outer_iterations",'i',&n_outer_iterations) == FETCH_SUCCESS,"reading number_of_outer_iterations from %s",input_file_path); print_continue("%i outer of %i inner iterations to be done",n_outer_iterations,n_inner_iterations); // read existing case and data case_file = fopen(case_file_path,"r"); if(case_file != NULL) { print_info("reading existing case file %s",case_file_path); read_case(case_file, &n_variables_old, &variable_order_old, &n_nodes, &node, &n_faces, &face, &n_elements, &element, &n_boundaries_old, &boundary_old); fclose(case_file); n = 0; for(i = 0; i < n_variables; i ++) n += n_elements*ORDER_TO_N_BASIS(variable_order[i]); if(fetch_value(input_file,"initial_data_file_path",'s',data_file_path) == FETCH_SUCCESS) { print_info("reading existing data file %s",data_file_path); exit_if_false(data_file = fopen(data_file_path,"r"),"opening %s",data_file_path); read_data(data_file, &n_u_old, &u_old, &outer_iteration); fclose(data_file); exit_if_false(n_u_old == n,"case and initial data does not match"); } } // construct new case else { print_info("reading the geometry file path"); exit_if_false(fetch_value(input_file,"geometry_file_path",'s',geometry_file_path) == FETCH_SUCCESS,"reading geometry_file_path from %s",input_file_path); print_continue(geometry_file_path); print_info("reading the geometry file %s",geometry_file_path); exit_if_false(geometry_file = fopen(geometry_file_path,"r"),"opening %s",geometry_file_path); read_geometry(geometry_file, &n_nodes, &node, &n_faces, &face, &n_elements, &element); fclose(geometry_file); print_continue("%i nodes, %i faces and %i elements",n_nodes,n_faces,n_elements); print_info("generating additional connectivity and geometry"); process_geometry(n_nodes, node, n_faces, face, n_elements, element); } // read the data file path and output frequency print_info("reading the data file path and output frequency"); exit_if_false(fetch_value(input_file,"data_file_path",'s',data_file_path) == FETCH_SUCCESS,"reading data_file_path from %s",input_file_path); if(fetch_value(input_file,"data_number_of_outer_iterations",'i',&data_n_outer_iterations) != FETCH_SUCCESS) data_n_outer_iterations = n_outer_iterations + outer_iteration; print_continue("data to be written to %s every %i outer iterations",data_file_path,data_n_outer_iterations); // read boundaries print_info("reading boundaries"); boundaries_input(input_file, n_faces, face, &n_boundaries, &boundary); print_continue("%i boundaries",n_boundaries); // read terms print_info("reading PDE terms"); terms_input(input_file, &n_terms, &term); print_continue("%i terms",n_terms); // update unknown indices and values print_info("updating the numbering of the degrees of freedom"); update_element_unknowns(n_variables_old, n_variables, variable_order_old, variable_order, n_elements, element, n_u_old, &n_u, u_old, &u); print_continue("%i degrees of freedom",n_u); // update face boundaries print_info("updating the face boundary associations"); update_face_boundaries(n_variables, n_faces, face, n_boundaries, boundary); // update integration print_info("updating integration"); i = update_face_integration(n_variables_old, n_variables, variable_order_old, variable_order, n_faces, face); if(i) print_continue("updated %i face quadratures",i); i = update_element_integration(n_variables_old, n_variables, variable_order_old, variable_order, n_elements, element); if(i) print_continue("updated %i element quadratures",i); // update numerics print_info("updating numerics"); i = update_face_numerics(n_variables_old, n_variables, variable_order_old, variable_order, n_faces, face, n_boundaries_old, boundary_old); if(i) print_continue("updated %i face interpolations",i); i = update_element_numerics(n_variables_old, n_variables, variable_order_old, variable_order, n_elements, element); if(i) print_continue("updated %i element interpolations",i); // write case file print_info("writing case file %s",case_file_path); exit_if_false(case_file = fopen(case_file_path,"w"),"opening %s",case_file_path); write_case(case_file, n_variables, variable_order, n_nodes, node, n_faces, face, n_elements, element, n_boundaries, boundary); fclose(case_file); // read the display file path and output frequency print_info("reading the display file path and output frequency"); if( fetch_value(input_file,"display_file_path",'s',display_file_path) == FETCH_SUCCESS && fetch_value(input_file,"display_number_of_outer_iterations",'i',&display_n_outer_iterations) == FETCH_SUCCESS ) print_continue("display to be written to %s every %i outer iterations",display_file_path,display_n_outer_iterations); else { display_n_outer_iterations = 0; warn_if_false(0,"display files will not be written"); } // initialise if(initial_input(input_file, n_variables, &initial)) { print_info("initialising the degrees of freedom"); initialise_values(n_variables, variable_order, n_elements, element, initial, u); } //-------------------------------------------------------------------// // allocate and initialise the system print_info("allocating and initialising the system"); SPARSE system = NULL; initialise_system(n_variables, variable_order, n_elements, element, n_u, &system); double *residual, *max_residual, *du, *max_du; exit_if_false(residual = (double *)malloc(n_u * sizeof(double)),"allocating the residuals"); exit_if_false(max_residual = (double *)malloc(n_variables * sizeof(double)),"allocating the maximum residuals"); exit_if_false(du = (double *)malloc(n_u * sizeof(double)),"allocating du"); exit_if_false(max_du = (double *)malloc(n_variables * sizeof(double)),"allocating the maximum changes"); exit_if_false(u_old = (double *)realloc(u_old, n_u * sizeof(double)),"re-allocating u_old"); timer_reset(); // iterate print_info("iterating"); n_outer_iterations += outer_iteration; for(; outer_iteration < n_outer_iterations; outer_iteration ++) { print_output("iteration %i", outer_iteration); for(i = 0; i < n_u; i ++) u_old[i] = u[i]; for(inner_iteration = 0; inner_iteration < n_inner_iterations; inner_iteration ++) { calculate_system(n_variables, variable_order, n_faces, face, n_elements, element, n_terms, term, n_u, u_old, u, system, residual); exit_if_false(sparse_solve(system, du, residual) == SPARSE_SUCCESS,"solving system"); for(i = 0; i < n_u; i ++) u[i] -= du[i]; calculate_maximum_changes_and_residuals(n_variables, variable_order, n_elements, element, du, max_du, residual, max_residual); for(i = 0; i < n_variables; i ++) sprintf(&print[26*i],"%.6e|%.6e ",max_du[i],max_residual[i]); print_continue("%s",print); } slope_limit(n_variables, variable_order, n_nodes, node, n_elements, element, n_boundaries, boundary, u); if(data_n_outer_iterations && outer_iteration % data_n_outer_iterations == 0) { generate_numbered_file_path(data_numbered_file_path, data_file_path, outer_iteration); print_info("writing data to %s",data_numbered_file_path); exit_if_false(data_file = fopen(data_numbered_file_path,"w"),"opening %s",data_numbered_file_path); write_data(data_file, n_u, u, outer_iteration); fclose(data_file); } if(display_n_outer_iterations && outer_iteration % display_n_outer_iterations == 0) { generate_numbered_file_path(display_numbered_file_path, display_file_path, outer_iteration); print_info("writing display to %s",data_numbered_file_path); exit_if_false(display_file = fopen(display_numbered_file_path,"w"),"opening %s",display_numbered_file_path); write_display(display_file, n_variables, variable_name, variable_order, n_elements, element, n_u, u); fclose(display_file); } timer_print(); } //-------------------------------------------------------------------// print_info("freeing all memory"); fclose(input_file); free(geometry_file_path); free(case_file_path); free(data_file_path); free(data_numbered_file_path); free(display_file_path); free(display_numbered_file_path); free(print); destroy_nodes(n_nodes,node); destroy_faces(n_faces,face,n_variables); destroy_elements(n_elements,element,n_variables); destroy_boundaries(n_boundaries_old,boundary_old); destroy_boundaries(n_boundaries,boundary); destroy_terms(n_terms,term); destroy_initial(n_variables,initial); free(variable_order_old); free(variable_order); destroy_matrix((void *)variable_name); free(u_old); free(u); sparse_destroy(system); free(residual); free(max_residual); free(du); free(max_du); return 0; }
void node_case_write(FILE *file, struct NODE *node) { exit_if_false(fwrite(node->x, sizeof(double), 2, file) == 2, "writing the node location"); }
void divergences_input(char *filename, char *constant, int *n_divergences, struct DIVERGENCE **divergence) { //open the file FILE *file = fopen(filename,"r"); exit_if_false(file != NULL,"opening input file"); //fetch the data FETCH fetch = fetch_new(DIVERGENCE_FORMAT,MAX_DIVERGENCES); exit_if_false(fetch != NULL,"allocating fetch"); int n_fetch = fetch_read(file,DIVERGENCE_LABEL,fetch); exit_if_false(n_fetch > 1,"no divergences found in input file"); warn_if_false(n_fetch < MAX_DIVERGENCES,"maximum number of divergences reached"); //allocate pointers struct DIVERGENCE *d = divergences_new(n_fetch,NULL); exit_if_false(d != NULL,"allocating divergences"); //counters int i, j, n = 0, info; //temporary storage char direction; int var_offset, dif_offset, pow_offset, dif[2]; char *var_string = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); char *dif_string = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); char *pow_string = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); char *temp = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); exit_if_false(var_string != NULL && dif_string != NULL && pow_string != NULL && temp != NULL,"allocating temporary strings"); int *vars = (int *)malloc(MAX_DIVERGENCE_VARIABLES * sizeof(int)); int *difs = (int *)malloc(MAX_DIVERGENCE_VARIABLES * sizeof(int)); int *pows = (int *)malloc(MAX_DIVERGENCE_VARIABLES * sizeof(int)); exit_if_false(vars != NULL && difs != NULL && pows != NULL,"allocating temporary data"); for(i = 0; i < n_fetch; i ++) { //equation fetch_get(fetch, i, 0, &d[n].equation); //direction fetch_get(fetch, i, 1, &direction); if(direction == 'x') { d[n].direction = 0; } else if(direction == 'y') { d[n].direction = 1; } else { warn_if_false(info = 0,"skipping divergence with unrecognised direction"); continue; } //get the variable, differential and power strings fetch_get(fetch, i, 2, var_string); fetch_get(fetch, i, 3, dif_string); fetch_get(fetch, i, 4, pow_string); for(j = 0; j < strlen(var_string); j ++) if(var_string[j] == ',') var_string[j] = ' '; for(j = 0; j < strlen(dif_string); j ++) if(dif_string[j] == ',') dif_string[j] = ' '; for(j = 0; j < strlen(pow_string); j ++) if(pow_string[j] == ',') pow_string[j] = ' '; //read each variable in turn var_offset = dif_offset = pow_offset = d[n].n_variables = 0; while(var_offset < strlen(var_string)) { info = 1; //read the variable index from the string info *= sscanf(&var_string[var_offset],"%s",temp) == 1; info *= sscanf(temp,"%i",&vars[d[n].n_variables]) == 1; var_offset += strlen(temp) + 1; //read the x and y differentials and convert to a differential index info *= sscanf(&dif_string[dif_offset],"%s",temp) == 1; j = dif[0] = dif[1] = 0; if(info) { while(temp[j] != '\0') { dif[0] += (temp[j] == 'x'); dif[1] += (temp[j] == 'y'); j ++; } difs[d[n].n_variables] = differential_index[dif[0]][dif[1]]; } dif_offset += strlen(temp) + 1; //read the variable powers from the string info *= sscanf(&pow_string[pow_offset],"%s",temp) == 1; info *= sscanf(temp,"%i",&pows[d[n].n_variables]) == 1; pow_offset += strlen(temp) + 1; warn_if_false(info,"skipping divergence with unrecognised variable format"); if(!info) continue; //next variable d[n].n_variables ++; } //allocate the variable and differential arrays d[n].variable = (int *)malloc(d[n].n_variables * sizeof(int)); //?? SHOULD BE MEMORY FUNCTIONS DOING THIS d[n].differential = (int *)malloc(d[n].n_variables * sizeof(int)); d[n].power = (int *)malloc(d[n].n_variables * sizeof(int)); exit_if_false(d[n].variable != NULL && d[n].differential != NULL && d[n].power != NULL,"allocating divergence variables and differentials"); //copy over for(j = 0; j < d[n].n_variables; j ++) { d[n].variable[j] = vars[j]; d[n].differential[j] = difs[j]; d[n].power[j] = pows[j]; } //implicit fraction fetch_get(fetch, i, 5, &d[n].implicit); //constant expression strcpy(temp, constant); j = strlen(temp); fetch_get(fetch, i, 6, &temp[j]); d[n].constant = expression_generate(temp); warn_if_false(d[n].constant != NULL,"skipping divergence for which the expression generation failed"); if(d[n].constant == NULL) continue; //printf("divergence %i expression -> ",n); expression_print(d[n].constant); printf("\n"); //increment the number of divergences n ++; } //check numbers fetch_destroy(fetch); fetch = fetch_new("",MAX_DIVERGENCES); warn_if_false(fetch_read(file,DIVERGENCE_LABEL,fetch) == n,"skipping divergences with unrecognised formats"); //copy over *n_divergences = n; *divergence = d; //clean up fclose(file); fetch_destroy(fetch); free(var_string); free(dif_string); free(pow_string); free(temp); free(vars); free(difs); free(pows); }
void zones_input(char *filename, int n_faces, struct FACE *face, int n_cells, struct CELL *cell, int *n_zones, struct ZONE **zone) { //counters int i, j, n = 0, info; //open the file FILE *file = fopen(filename,"r"); exit_if_false(file != NULL,"opening input file"); //fetch the data from the file FETCH fetch = fetch_new(ZONE_FORMAT, MAX_ZONES); exit_if_false(fetch != NULL,"allocating zone input"); int n_fetch = fetch_read(file, ZONE_LABEL, fetch); exit_if_false(n_fetch > 1,"no zones found in input file"); warn_if_false(n_fetch < MAX_ZONES,"maximum number of zones reached"); //allocate zones struct ZONE *z = zones_new(n_fetch, NULL); exit_if_false(z != NULL,"allocating zones"); //temporary storage int offset, index[2]; char *range = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); char *temp = (char *)malloc(MAX_STRING_LENGTH * sizeof(char)); exit_if_false(range != NULL && temp != NULL,"allocating temporary storage"); //consider each feteched line for(i = 0; i < n_fetch; i ++) { //get zone parameters fetch_get(fetch, i, 0, &z[n].location); fetch_get(fetch, i, 2, &z[n].variable); fetch_get(fetch, i, 3, z[n].condition); fetch_get(fetch, i, 4, &z[n].value); //get the range string fetch_get(fetch, i, 1, range); //convert comma delimiters to whitespace for(j = 0; j < strlen(range); j ++) if(range[j] == ',') range[j] = ' '; //sequentially read ranges offset = info = 0; while(offset < strlen(range)) { //read the range from the string info = sscanf(&range[offset],"%s",temp) == 1; info *= sscanf(temp,"%i:%i",&index[0],&index[1]) == 2; warn_if_false(info,"skipping zone with unrecognised range"); if(!info) break; //store zone in the elements in the range if(z[n].location == 'f') for(j = index[0]; j <= index[1]; j ++) exit_if_false(face_zone_add(&face[j],&z[n]),"adding a face zone"); if(z[n].location == 'c') for(j = index[0]; j <= index[1]; j ++) exit_if_false(cell_zone_add(&cell[j],&z[n]),"adding a cell zone"); //move to the next range in the string offset += strlen(temp) + 1; } //increment zone n += info; } //check numbers fetch_destroy(fetch); fetch = fetch_new("",MAX_ZONES); warn_if_false(fetch_read(file,ZONE_LABEL,fetch) == n,"skipping zones with unrecognised formats"); //resize zone list struct ZONE *z_new = zones_new(n, z); exit_if_false(zone != NULL,"re-allocating zones"); for(i = 0; i < n_faces; i ++) for(j = 0; j < face[i].n_zones; j ++) face[i].zone[j] += z_new - z; for(i = 0; i < n_cells; i ++) for(j = 0; j < cell[i].n_zones; j ++) cell[i].zone[j] += z_new - z; //copy over *n_zones = n; *zone = z_new; //clean up fclose(file); fetch_destroy(fetch); free(range); free(temp); }
void node_case_get(FILE *file, struct NODE *node) { exit_if_false(fread(node->x, sizeof(double), 2, file) == 2, "reading the node location"); }
void read_geometry(char *filename, int *n_nodes, struct NODE **node, int *n_faces, struct FACE **face, int *n_cells, struct CELL **cell) { FILE *file = fopen(filename,"r"); exit_if_false(file != NULL,"opening geometry file"); char *line; exit_if_false(allocate_character_vector(&line, MAX_STRING_LENGTH),"allocating line string"); int i; while(fgets(line, MAX_STRING_LENGTH, file) != NULL) { if(strncmp(line,"NODES",5) == 0) { exit_if_false(sscanf(&line[6],"%i",n_nodes) == 1,"reading the number of nodes"); *node = nodes_new(*n_nodes, NULL); exit_if_false(*node != NULL,"allocating the nodes"); for(i = 0; i < *n_nodes; i ++) node_geometry_get(file, &(*node)[i]); } if(strncmp(line,"FACES",5) == 0) { exit_if_false(sscanf(&line[6],"%i",n_faces) == 1,"reading the number of faces"); *face = faces_new(*n_faces, NULL); exit_if_false(*face != NULL,"allocating the faces"); for(i = 0; i < *n_faces; i ++) face_geometry_get(file, &(*face)[i], *node); } if(strncmp(line,"CELLS",5) == 0) { exit_if_false(sscanf(&line[6],"%i",n_cells) == 1,"reading the number of cells"); *cell = cells_new(*n_cells, NULL); exit_if_false(*cell != NULL,"allocating the cells"); for(i = 0; i < *n_cells; i ++) cell_geometry_get(file, &(*cell)[i], *face); } } exit_if_false(*n_nodes > 0,"finding nodes in the geometry file"); exit_if_false(*n_faces > 0,"finding faces in the geometry file"); exit_if_false(*n_cells > 0,"finding cells in the geometry file"); free(line); fclose(file); }
void node_geometry_get(FILE *file, struct NODE *node) { int info = fscanf(file,"%lf %lf\n",&(node->x[0]),&(node->x[1])); exit_if_false(info == 2, "reading a node's coordinates"); }
void cell_case_get(FILE *file, int n_variables, struct FACE *face, struct CELL *cell, struct ZONE *zone) { int i, n, *index; index = (int *)malloc(MAX_INDICES * sizeof(int)); exit_if_false(index != NULL,"allocating temporary storage"); exit_if_false(fread(&(cell->n_faces), sizeof(int), 1, file) == 1, "reading the number of cell faces"); exit_if_false(cell_face_new(cell),"allocating cell faces"); exit_if_false(fread(index, sizeof(int), cell->n_faces, file) == cell->n_faces, "reading the cell faces"); for(i = 0; i < cell->n_faces; i ++) cell->face[i] = &face[index[i]]; exit_if_false(cell_oriented_new(cell),"allocating cell orientations"); exit_if_false(fread(cell->oriented, sizeof(int), cell->n_faces, file) == cell->n_faces, "reading cell orientations"); exit_if_false(fread(&(cell->area), sizeof(double), 1, file) == 1, "reading the cell area"); exit_if_false(fread(cell->centroid, sizeof(double), 2, file) == 2, "reading the cell centroid"); exit_if_false(fread(&(cell->n_zones), sizeof(int), 1, file) == 1, "reading the number of cell zones"); exit_if_false(cell_zone_new(cell),"allocating cell zones"); exit_if_false(fread(index, sizeof(int), cell->n_zones, file) == cell->n_zones, "reading the cell zones"); for(i = 0; i < cell->n_zones; i ++) cell->zone[i] = &zone[index[i]]; exit_if_false(cell_order_new(n_variables,cell),"allocating cell orders"); exit_if_false(fread(cell->order, sizeof(int), n_variables, file) == n_variables, "reading the cell orders"); exit_if_false(cell_n_stencil_new(n_variables,cell),"allocating cell stencil sizes"); exit_if_false(fread(cell->n_stencil, sizeof(int), n_variables, file) == n_variables, "reading the cell stencil sizes"); exit_if_false(cell_stencil_new(n_variables,cell),"allocating cell stencils"); exit_if_false(cell_matrix_new(n_variables,cell),"allocating cell matrices"); for(i = 0; i < n_variables; i ++) { exit_if_false(fread(cell->stencil[i], sizeof(int), cell->n_stencil[i], file) == cell->n_stencil[i],"reading the cell stencil"); n = ORDER_TO_POWERS(cell->order[i]) * cell->n_stencil[i]; exit_if_false(fread(cell->matrix[i][0], sizeof(double), n, file) == n,"reading the cell matrix"); } free(index); }
void calculate_cell_reconstruction_matrices(int n_variables, double *weight_exponent, int *maximum_order, struct FACE *face, int n_cells, struct CELL *cell, struct ZONE *zone) { int c, u, i, j, k, l; int order, n_powers, n_stencil; //find the overall maximum order int maximum_maximum_order = 0; for(u = 0; u < n_variables; u ++) if(maximum_order[u] > maximum_maximum_order) maximum_maximum_order = maximum_order[u]; //cell structure allocation for(c = 0; c < n_cells; c ++) exit_if_false(cell_matrix_new(n_variables, &cell[c]),"allocating cell matrices"); //numerics values double **matrix, *weight; int n_constraints, *constraint; exit_if_false(allocate_double_matrix(&matrix,ORDER_TO_POWERS(maximum_maximum_order),MAX_STENCIL),"allocating matrix"); exit_if_false(allocate_integer_vector(&constraint,MAX_STENCIL),"allocating constraints"); exit_if_false(allocate_double_vector(&weight,MAX_STENCIL),"allocating weights"); //stencil element properties int s_id, s_index; struct ZONE *s_zone; char s_location, *s_condition; double s_area, *s_centroid, s_weight; //integration double x[2]; int differential[2], d; //CV polygon int n_polygon; double ***polygon; exit_if_false(allocate_double_pointer_matrix(&polygon,MAX(MAX_CELL_FACES,4),2),"allocating polygon memory"); for(c = 0; c < n_cells; c ++) { for(u = 0; u < n_variables; u ++) { //problem size order = cell[c].order[u]; n_powers = ORDER_TO_POWERS(order); n_stencil = cell[c].n_stencil[u]; n_constraints = 0; for(i = 0; i < n_stencil; i ++) { //stencil element properties s_id = cell[c].stencil[u][i]; s_index = ID_TO_INDEX(s_id); s_zone = &zone[ID_TO_ZONE(s_id)]; s_location = s_zone->location; s_condition = s_zone->condition; if(s_location == 'f') { s_centroid = face[s_index].centroid; s_area = face[s_index].area; } else if(s_location == 'c') { s_centroid = cell[s_index].centroid; s_area = cell[s_index].area; } else exit_if_false(0,"recognising zone location"); s_weight = (s_centroid[0] - cell[c].centroid[0])*(s_centroid[0] - cell[c].centroid[0]); s_weight += (s_centroid[1] - cell[c].centroid[1])*(s_centroid[1] - cell[c].centroid[1]); s_weight = 1.0/pow(s_weight,0.5*weight_exponent[u]); if(s_location == 'c' && s_index == c) s_weight = 1.0; weight[i] = s_weight; //unknown and dirichlet conditions have zero differentiation differential[0] = differential[1] = 0; //other conditions have differential determined from numbers of x and y-s in the condition string if(s_condition[0] != 'u' && s_condition[0] != 'd') { j = 0; while(s_condition[j] != '\0') { differential[0] += (s_condition[j] == 'x'); differential[1] += (s_condition[j] == 'y'); j ++; } } //index for the determined differential d = differential_index[differential[0]][differential[1]]; //unknowns if(s_condition[0] == 'u') { /*//fit unknowns to centroid points x[0] = s_centroid[0] - cell[c].centroid[0]; x[1] = s_centroid[1] - cell[c].centroid[1]; for(j = 0; j < n_powers; j ++) { matrix[j][i] = polynomial_coefficient[d][j]* integer_power(x[0],polynomial_power_x[d][j])* integer_power(x[1],polynomial_power_y[d][j])* s_weight; }*/ //fit unknowns to CV average n_polygon = generate_control_volume_polygon(polygon, s_index, s_location, face, cell); for(j = 0; j < n_powers; j ++) matrix[j][i] = 0.0; for(j = 0; j <= order; j ++) { for(k = 0; k < n_polygon; k ++) { x[0] = 0.5*polygon[k][0][0]*(1.0 - gauss_x[order][j]) + 0.5*polygon[k][1][0]*(1.0 + gauss_x[order][j]) - cell[c].centroid[0]; x[1] = 0.5*polygon[k][0][1]*(1.0 - gauss_x[order][j]) + 0.5*polygon[k][1][1]*(1.0 + gauss_x[order][j]) - cell[c].centroid[1]; for(l = 0; l < n_powers; l ++) { //[face integral of polynomial integrated wrt x] * [x normal] / [CV area] matrix[l][i] += polynomial_coefficient[d][l] * (1.0 / (polynomial_power_x[d][l] + 1.0)) * integer_power(x[0],polynomial_power_x[d][l]+1) * integer_power(x[1],polynomial_power_y[d][l]) * s_weight * gauss_w[order][j] * 0.5 * (polygon[k][1][1] - polygon[k][0][1]) / s_area; } } } } //knowns else { //known faces fit to face average if(s_location == 'f') { for(j = 0; j < n_powers; j ++) matrix[j][i] = 0.0; for(j = 0; j < order; j ++) { x[0] = 0.5*face[s_index].node[0]->x[0]*(1.0 - gauss_x[order-1][j]) + 0.5*face[s_index].node[1]->x[0]*(1.0 + gauss_x[order-1][j]) - cell[c].centroid[0]; x[1] = 0.5*face[s_index].node[0]->x[1]*(1.0 - gauss_x[order-1][j]) + 0.5*face[s_index].node[1]->x[1]*(1.0 + gauss_x[order-1][j]) - cell[c].centroid[1]; for(k = 0; k < n_powers; k ++) { matrix[k][i] += polynomial_coefficient[d][k] * integer_power(x[0],polynomial_power_x[d][k]) * integer_power(x[1],polynomial_power_y[d][k]) * s_weight*gauss_w[order-1][j]*0.5; } } } //cells need implementing //if(s_location == 'c') //{ //} } //constraints are the centre cell and any dirichlet boundaries if((s_location == 'c' && s_index == c) || s_condition[0] == 'd') constraint[n_constraints++] = i; } //solve if(n_constraints > 0) exit_if_false(constrained_least_squares(n_stencil,n_powers,matrix,n_constraints,constraint) == LS_SUCCESS, "doing CLS calculation"); else exit_if_false(least_squares(n_stencil,n_powers,matrix) == LS_SUCCESS,"doing LS calculation"); //multiply by the weights for(i = 0; i < n_powers; i ++) for(j = 0; j < n_stencil; j ++) matrix[i][j] *= weight[j]; //store in the cell structure for(i = 0; i < n_powers; i ++) for(j = 0; j < n_stencil; j ++) cell[c].matrix[u][i][j] = matrix[i][j]; } } //clean up free_matrix((void**)matrix); free_vector(constraint); free_vector(weight); free_matrix((void**)polygon); }