static void gfs_refine_surface_refine (GfsRefine * refine, GfsSimulation * sim) { RefineCut p; p.refine = refine; p.domain = GFS_DOMAIN (sim); p.surface = GFS_REFINE_SURFACE (refine)->surface; if (GFS_SURFACE (p.surface)->twod) { if (GFS_SURFACE (p.surface)->s) gfs_domain_traverse_cut_2D (GFS_DOMAIN (sim), GFS_REFINE_SURFACE (refine)->surface, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, (FttCellTraverseCutFunc) refine_cut_cell, &p); else g_assert_not_implemented (); } else { if (GFS_SURFACE (p.surface)->s) gfs_domain_traverse_cut (GFS_DOMAIN (sim), GFS_REFINE_SURFACE (refine)->surface, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, (FttCellTraverseCutFunc) refine_cut_cell, &p); else gfs_domain_cell_traverse (GFS_DOMAIN (sim), FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) refine_implicit_cell, &p); } }
static void wave_init (GfsWave * wave) { wave->nk = 25; wave->ntheta = 24; wave->alpha_s = 0.; /* default for g is acceleration of gravity on Earth with kilometres as spatial units, hours as time units and Hz as frequency units */ GFS_SIMULATION (wave)->physical_params.g = 9.81/1000.*3600.; GfsAdvectionParams * par = &GFS_SIMULATION (wave)->advection_params; par->gradient = gfs_center_van_leer_gradient; par->flux = gfs_face_advection_flux; par->use_centered_velocity = FALSE; static GfsDerivedVariableInfo derived_variable[] = { { "Hs", "Significant wave height", cell_hs }, { "Energy", "Wave energy", cell_E }, { "Frequency", "Wave frequency", cell_frequency }, { "Direction", "Wave direction (angle)", cell_direction }, { NULL, NULL, NULL} }; GfsDerivedVariableInfo * v = derived_variable; while (v->name) { g_assert (gfs_domain_add_derived_variable (GFS_DOMAIN (wave), *v)); v++; } }
static void wave_read (GtsObject ** o, GtsFile * fp) { (* GTS_OBJECT_CLASS (gfs_wave_class ())->parent_class->read) (o, fp); if (fp->type == GTS_ERROR) return; GfsWave * wave = GFS_WAVE (*o); if (fp->type == '{') { GtsFileVariable var[] = { {GTS_UINT, "nk", TRUE, &wave->nk}, {GTS_UINT, "ntheta", TRUE, &wave->ntheta}, {GTS_DOUBLE, "alpha_s", TRUE, &wave->alpha_s}, {GTS_NONE} }; gts_file_assign_variables (fp, var); if (fp->type == GTS_ERROR) return; } GfsDomain * domain = GFS_DOMAIN (wave); guint ik, ith; wave->F = gfs_matrix_new (wave->nk, wave->ntheta, sizeof (GfsVariable *)); for (ik = 0; ik < wave->nk; ik++) for (ith = 0; ith < wave->ntheta; ith++) { gchar * name = g_strdup_printf ("F%d_%d", ik, ith); gchar * description = g_strdup_printf ("Action density for f = %g Hz and theta = %g degrees", frequency (ik), theta (ith, wave->ntheta)*180./M_PI); wave->F[ik][ith] = gfs_domain_get_or_add_variable (domain, name, description); g_assert (wave->F[ik][ith]); g_free (name); g_free (description); } }
static gboolean gfs_init_wave_event (GfsEvent * event, GfsSimulation * sim) { if ((* GFS_EVENT_CLASS (GTS_OBJECT_CLASS (gfs_init_wave_class ())->parent_class)->event) (event, sim)) { gfs_catch_floating_point_exceptions (); gfs_domain_cell_traverse (GFS_DOMAIN (sim), FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) init_energy, event); gfs_restore_fpe_for_function (GFS_INIT_WAVE (event)->d); gfs_catch_floating_point_exceptions (); gfs_domain_cell_traverse (GFS_DOMAIN (sim), FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) scale_energy, event); gfs_restore_fpe_for_function (GFS_INIT_WAVE (event)->hs); return TRUE; } return FALSE; }
static void refine_solid_destroy (GtsObject * object) { gfs_domain_remove_derived_variable (GFS_DOMAIN (gfs_object_simulation (object)), "SolidCurvature"); (* GTS_OBJECT_CLASS (gfs_refine_solid_class ())->parent_class->destroy) (object); }
static void swap_face_fractions (GfsSimulation * sim) { GfsDomain * domain = GFS_DOMAIN (sim); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_ALL, -1, (FttCellTraverseFunc) swap_fractions, GFS_SIMULATION_MOVING (sim)->old_solid); gts_container_foreach (GTS_CONTAINER (domain), (GtsFunc) foreach_box, NULL); }
static gboolean output_spectra_event (GfsEvent * event, GfsSimulation * sim) { if ((* GFS_EVENT_CLASS (GTS_OBJECT_CLASS (gfs_output_spectra_class ())->parent_class)->event) (event, sim)) { GfsDomain * domain = GFS_DOMAIN (sim); GfsOutputSpectra * v = GFS_OUTPUT_SPECTRA (event); fftw_plan p; Datawrite data; data.fp = GFS_OUTPUT (event)->file->fp; data.L = v->L; data.kmax = init_kmax(v->L); data.dir1 = v->dir[0]; data.dir2 = v->dir[1]; fill_cartesian_matrix( v->cgd, v->v, domain); switch (v->Ndim) { case 1: { data.n1 = ( v->cgd->n[v->dir[0]] / 2 ) + 1; data.out = fftw_malloc( sizeof(fftw_complex)*data.n1 ); p = fftw_plan_dft_r2c_1d( v->cgd->n[v->dir[0]], v->cgd->v, data.out, FFTW_ESTIMATE); fftw_execute(p); write_spectra_1D ( &data ); break; } case 2: { data.n1 = v->cgd->n[v->dir[0]]; data.n2 = ( v->cgd->n[v->dir[1]] / 2 ) + 1; data.out = fftw_malloc( sizeof(fftw_complex)*v->cgd->n[v->dir[0]]*data.n2 ); p = fftw_plan_dft_r2c_2d( v->cgd->n[v->dir[0]], v->cgd->n[v->dir[1]], v->cgd->v, data.out, FFTW_ESTIMATE); fftw_execute(p); write_spectra_2D ( &data ); break; } case 3: { data.n1 = v->cgd->n[0]; data.n2 = v->cgd->n[1]; data.n3 = ( v->cgd->n[2] / 2 ) + 1; data.out = fftw_malloc( sizeof(fftw_complex)*v->cgd->n[0]*v->cgd->n[1]*data.n3 ); p = fftw_plan_dft_r2c_3d( v->cgd->n[0], v->cgd->n[1], v->cgd->n[2], v->cgd->v, data.out, FFTW_ESTIMATE); fftw_execute(p); write_spectra_3D ( &data ); break; } default: g_assert_not_reached (); } fftw_destroy_plan(p); fftw_free ( data.out ); return TRUE; } return FALSE; }
static void refine_distance_destroy (GtsObject * object) { GfsRefineDistance * d = GFS_REFINE_DISTANCE (object); if (d->stree) gts_bb_tree_destroy (d->stree, TRUE); gfs_domain_remove_derived_variable (GFS_DOMAIN (gfs_object_simulation (object)), "Distance"); (* GTS_OBJECT_CLASS (gfs_refine_distance_class ())->parent_class->destroy) (object); }
static gboolean gfs_init_stokes_wave_event (GfsEvent * event, GfsSimulation * sim) { if ((* GFS_EVENT_CLASS (GTS_OBJECT_CLASS (gfs_init_stokes_wave_class ())->parent_class)->event) (event, sim)) { GfsVariable ** velocity = gfs_domain_velocity (GFS_DOMAIN (sim)); GfsVariable * t = gfs_variable_from_name (GFS_DOMAIN (sim)->variables, "T"); g_assert (velocity); g_assert (t); gfs_domain_cell_traverse (GFS_DOMAIN (sim), FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) init_velocity, velocity); GfsSurface * surface = GFS_SURFACE (gts_object_new (GTS_OBJECT_CLASS (gfs_surface_class ()))); surface->f = gfs_function_spatial_new (gfs_function_spatial_class (), stokes_height); gfs_object_simulation_set (surface->f, sim); gfs_domain_init_fraction (GFS_DOMAIN (sim), GFS_GENERIC_SURFACE (surface), t); gts_object_destroy (GTS_OBJECT (surface)); return TRUE; } return FALSE; }
static void refine_solid_read (GtsObject ** o, GtsFile * fp) { GfsRefineSolid * refine = GFS_REFINE_SOLID (*o); GfsDerivedVariableInfo v = { "SolidCurvature", "curvature of the solid boundary", solid_curvature }; refine->v = gfs_domain_add_derived_variable (GFS_DOMAIN (gfs_object_simulation (*o)), v); if (!refine->v) { gts_file_error (fp, "derived variable `SolidCurvature' already defined"); return; } (* GTS_OBJECT_CLASS (gfs_refine_solid_class ())->parent_class->read) (o, fp); }
static void scale_energy (FttCell * cell, GfsInitWave * event) { GfsWave * wave = GFS_WAVE (gfs_object_simulation (event)); gdouble E = cell_E (cell, NULL, GFS_DOMAIN (wave)); if (E > 0.) { gdouble Hs = gfs_function_value (event->hs, cell); gdouble scaling = Hs*Hs/(16.*E); guint ik, ith; for (ik = 0; ik < wave->nk; ik++) for (ith = 0; ith < wave->ntheta; ith++) GFS_VALUE (cell, wave->F[ik][ith]) *= scaling; } }
static void gfs_refine_solid_refine (GfsRefine * refine, GfsSimulation * sim) { if (sim->solids) { RefineCut p; p.refine = refine; p.domain = GFS_DOMAIN (sim); GSList * i = sim->solids->items; while (i) { p.surface = GFS_SOLID (i->data)->s; gfs_catch_floating_point_exceptions (); if (GFS_SURFACE (p.surface)->s) gfs_domain_traverse_cut (GFS_DOMAIN (sim), p.surface, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, (FttCellTraverseCutFunc) refine_cut_cell, &p); else gfs_domain_cell_traverse (GFS_DOMAIN (sim), FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) refine_implicit_cell, &p); gfs_restore_fpe_for_function (refine->maxlevel); i = i->next; } } }
static void gfs_source_darcy_read (GtsObject ** o, GtsFile * fp) { (* GTS_OBJECT_CLASS (gfs_source_darcy_class ())->parent_class->read) (o, fp); if (fp->type == GTS_ERROR) return; printf("\ntesting1\n"); /*Check if source darcy has already been added or not*/ FttComponent c; for (c = 0; c < FTT_DIMENSION; c++) { GfsVariable * v = GFS_SOURCE_VELOCITY (*o)->v[c]; if (v->sources) { GSList * i = GTS_SLIST_CONTAINER (v->sources)->items; while (i) { if (i->data != *o && GFS_IS_SOURCE_DARCY (i->data)) { gts_file_error (fp, "variable '%s' cannot have multiple Darcy source terms", v->name); return; } i = i->next; } } } printf("\ntesting2\n"); GfsDomain * domain = GFS_DOMAIN (gfs_object_simulation (*o)); GfsSourceDarcy * s = GFS_SOURCE_DARCY (*o); printf("\ntesting3\n"); s->darcycoeff = gfs_function_new (gfs_function_class (), 0.); gfs_function_read (s->darcycoeff, gfs_object_simulation (s), fp); printf("\ntesting4\n"); if (fp->type != '\n') { s->forchhicoeff = gfs_function_new (gfs_function_class (), 0.); gfs_function_read (s->forchhicoeff, gfs_object_simulation (s), fp); } if (s->beta < 1.) for (c = 0; c < FTT_DIMENSION; c++) s->u[c] = gfs_temporary_variable (domain); else { GFS_SOURCE_GENERIC (s)->centered_value = NULL; GFS_SOURCE_GENERIC (s)->mac_value = NULL; } printf("\ntesting5\n"); }
static void gfs_init_wave_read (GtsObject ** o, GtsFile * fp) { (* GTS_OBJECT_CLASS (gfs_init_wave_class ())->parent_class->read) (o, fp); if (fp->type == GTS_ERROR) return; GfsDomain * domain = GFS_DOMAIN (gfs_object_simulation (*o)); if (!GFS_IS_WAVE (domain)) { gts_file_error (fp, "GfsInitWave can only be used within a GfsWave simulation"); return; } gfs_function_read (GFS_INIT_WAVE (*o)->d, domain, fp); if (fp->type == GTS_ERROR) return; gfs_function_read (GFS_INIT_WAVE (*o)->hs, domain, fp); }
static void gfs_skew_symmetric_momentum (GfsSimulation * sim, FaceData * fd, GfsVariable **gmac) { GfsDomain * domain = GFS_DOMAIN (sim); FttComponent c; FttDirection d; /* it is used for implementation of viscosity (improve?) */ GfsSourceDiffusion * dif = source_diffusion_viscosity (fd->u[0]); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) advance_face_values, fd); /* boundary conditions */ for (d = 0; d < FTT_NEIGHBORS; d++) gfs_domain_bc (domain, FTT_TRAVERSE_LEAFS, -1, fd->velfaces[d]); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) advection_term, fd); if (dif) { DataDif dd = { dif , sim->physical_params.alpha, fd }; gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) diffusion_term, &dd); } /* regularize flux at faces */ for (c = 0; c < FTT_DIMENSION; c++) gfs_domain_face_bc (domain, c, fd->u[c]); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) obtain_face_fluxes, NULL); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) update_vel, fd); gfs_velocity_face_sources (domain, fd->u, (*fd->dt), sim->physical_params.alpha, gmac); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) correct_face_velocity, NULL); }
static gboolean gfs_init_face_values_event (GfsEvent * event, GfsSimulation * sim) { if ((* GFS_EVENT_CLASS (GTS_OBJECT_CLASS (gfs_init_face_values_class ())->parent_class)->event) (event, sim)) { GSList * i = GFS_INIT (event)->f; while (i) { VarFunc * vf = i->data; FaceInitData data; FttComponent c = FTT_DIMENSION; if (!strcmp (vf->v[0]->name, "U")) { if (vf->n > 1) g_assert_not_implemented (); data.v1 = GFS_SKEW_SYMMETRIC(sim)->velfaces[0]; data.v2 = GFS_SKEW_SYMMETRIC(sim)->velfaces[1]; data.f = vf->f[0]; c = FTT_X; } else if (!strcmp (vf->v[0]->name, "V")) { data.v1 = GFS_SKEW_SYMMETRIC(sim)->velfaces[2]; data.v2 = GFS_SKEW_SYMMETRIC(sim)->velfaces[3]; data.f = vf->f[0]; c = FTT_Y; } #if (!FTT_2D) else if (!strcmp (vf->v[0]->name, "W")) { data.v1 = GFS_SKEW_SYMMETRIC(sim)->velfaces[4]; data.v2 = GFS_SKEW_SYMMETRIC(sim)->velfaces[5]; data.f = vf->f[0]; c = FTT_Z; } #endif if (c < FTT_DIMENSION) { gfs_catch_floating_point_exceptions (); gfs_domain_face_traverse (GFS_DOMAIN (sim), c, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttFaceTraverseFunc) init_fd, &data); gfs_restore_fpe_for_function (vf->f[0]); } i = i->next; } return TRUE; } return FALSE; }
static void refine_height_read (GtsObject ** o, GtsFile * fp) { GfsDerivedVariableInfo v = { "Height", "vertical distance to the surface", cell_height }; v.data = *o; if (!gfs_domain_add_derived_variable (GFS_DOMAIN (gfs_object_simulation (*o)), v)) { gts_file_error (fp, "derived variable `Height' already defined"); return; } (* GTS_OBJECT_CLASS (gfs_refine_height_class ())->parent_class->read) (o, fp); if (fp->type == GTS_ERROR) return; if (!GFS_SURFACE (GFS_REFINE_SURFACE (*o)->surface)->s) { gts_file_error (fp, "RefineHeight only works with GTS surfaces"); return; } }
static gboolean gfs_source_darcy_event (GfsEvent * event, GfsSimulation * sim) { if ((* GFS_EVENT_CLASS (GTS_OBJECT_CLASS (gfs_source_darcy_class ())->parent_class)->event) (event, sim)) { /* do object-specific event here */ if (GFS_SOURCE_DARCY (event)->beta < 1.) { gfs_catch_floating_point_exceptions (); gfs_domain_traverse_layers (GFS_DOMAIN (sim), (FttCellTraverseFunc) save_darcy, event); if (gfs_restore_floating_point_exceptions ()) { GfsSourceDarcy * c = GFS_SOURCE_DARCY (event); gchar * s = g_strconcat ("\n", gfs_function_description (c->darcycoeff, FALSE), NULL); if (c->forchhicoeff) s = g_strconcat (s, "\n", gfs_function_description (c->forchhicoeff, FALSE), NULL); /* fixme: memory leaks */ g_message ("floating-point exception in user-defined function(s):%s", s); exit (1); } } return TRUE; } return FALSE; }
static void refine_distance_read (GtsObject ** o, GtsFile * fp) { GfsDerivedVariableInfo v = { "Distance", "distance to the surface", cell_distance }; v.data = *o; if (!gfs_domain_add_derived_variable (GFS_DOMAIN (gfs_object_simulation (*o)), v)) { gts_file_error (fp, "derived variable `Distance' already defined"); return; } (* GTS_OBJECT_CLASS (gfs_refine_distance_class ())->parent_class->read) (o, fp); if (fp->type == GTS_ERROR) return; GtsSurface * s = GFS_SURFACE (GFS_REFINE_SURFACE (*o)->surface)->s; if (!s) { gts_file_error (fp, "RefineDistance only works with GTS surfaces"); return; } GFS_REFINE_DISTANCE (*o)->stree = gts_bb_tree_surface (s); }
static void gfs_skew_symmetric_init (GfsSkewSymmetric * object) { object->beta = 0.05; GfsDomain * domain = GFS_DOMAIN (object); FttDirection d; for (d = 0; d < FTT_NEIGHBORS; d++) { gchar * name = g_strdup_printf ("Uface%d", d); gchar * descr = g_strdup_printf ("%d-component of face velocity", d); object->velfaces[d] = gfs_domain_add_variable (domain, name, descr); object->velfaces[d]->units = 1.; g_free(name); g_free(descr); name = g_strdup_printf ("Ufaceold%d", d); descr = g_strdup_printf ("%d-component of old face velocity", d); object->velold[d] = gfs_domain_add_variable (domain, name, descr); object->velold[d]->units = 1.; g_free(name); g_free(descr); } // gfs_variable_set_vector (object->velfaces, FTT_NEIGHBORS); // gfs_variable_set_vector (object->velold, FTT_NEIGHBORS); }
static void darcy_coefficients (GfsSourceDarcy * sd, FttCell * cell, GfsVariable ** u, FttComponent c1, gdouble f[2]) { GfsPorous * por = GFS_POROUS (gfs_object_simulation (sd)); GfsDomain * domain = GFS_DOMAIN(por); GfsFunction * porosity = por->porosity; GfsFunction * K = por->K; gdouble viscosity = 0.001; GfsVariable ** U = gfs_domain_velocity (domain); GfsSourceVelocity * sv = GFS_SOURCE_VELOCITY (domain); GfsSourceDiffusion * d = source_diffusion_viscosity (U[0]); if (d) viscosity = gfs_diffusion_cell (d->D, cell); f[0] = (viscosity * gfs_function_value (porosity,cell))/gfs_function_value (K,cell); gdouble modu; modu = fabs(sqrt(GFS_VALUE (cell, sv->v[0])*GFS_VALUE (cell, sv->v[0]) + GFS_VALUE (cell, sv->v[1])*GFS_VALUE (cell, sv->v[1]))); f[1] = (gfs_function_value (porosity,cell))/sqrt(gfs_function_value (K,cell))*modu; }
static void swap_face_fractions_back (GfsSimulation * sim) { gfs_domain_cell_traverse (GFS_DOMAIN (sim), FTT_PRE_ORDER, FTT_TRAVERSE_ALL, -1, (FttCellTraverseFunc) swap_fractions_back, GFS_SIMULATION_MOVING (sim)->old_solid); }
static void gfs_porous_run (GfsSimulation * sim) { GfsVariable * p, * pmac, * res = NULL, * g[FTT_DIMENSION], * gmac[FTT_DIMENSION]; GfsVariable ** gc = sim->advection_params.gc ? g : NULL; GfsDomain * domain; GfsPorous *por; GSList * i; domain = GFS_DOMAIN (sim); por = GFS_POROUS (sim); p = gfs_variable_from_name (domain->variables, "P"); g_assert (p); pmac = gfs_variable_from_name (domain->variables, "Pmac"); g_assert (pmac); FttComponent c; for (c = 0; c < FTT_DIMENSION; c++) { gmac[c] = gfs_temporary_variable (domain); if (sim->advection_params.gc) g[c] = gfs_temporary_variable (domain); else g[c] = gmac[c]; } gfs_variable_set_vector (gmac, FTT_DIMENSION); gfs_variable_set_vector (g, FTT_DIMENSION); gfs_simulation_refine (sim); gfs_simulation_init (sim); i = domain->variables; while (i) { if (GFS_IS_VARIABLE_RESIDUAL (i->data)) res = i->data; i = i->next; } gfs_simulation_set_timestep (sim); if (sim->time.i == 0) { /*inserted changes inside this function*/ gfs_approximate_projection_por (domain, por, &sim->approx_projection_params, sim->advection_params.dt, p, sim->physical_params.alpha, res, g, NULL); gfs_simulation_set_timestep (sim); gfs_advance_tracers (sim, sim->advection_params.dt/2.); } else if (sim->advection_params.gc) gfs_update_gradients_por (domain, por, p, sim->physical_params.alpha, g); while (sim->time.t < sim->time.end && sim->time.i < sim->time.iend) { gdouble tstart = gfs_clock_elapsed (domain->timer); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_do, sim); /*inserted changes */ gfs_pre_projection (domain, por, FTT_DIMENSION); if (sim->advection_params.linear) { /* linearised advection */ gfs_domain_face_traverse (domain, FTT_XYZ, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttFaceTraverseFunc) gfs_face_reset_normal_velocity, NULL); gfs_domain_face_traverse (domain, FTT_XYZ, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttFaceTraverseFunc) gfs_face_interpolated_normal_velocity, sim->u0); } else gfs_predicted_face_velocities (domain, FTT_DIMENSION, &sim->advection_params); gfs_variables_swap (p, pmac); gfs_mac_projection_por (domain, por, &sim->projection_params, sim->advection_params.dt/2., p, sim->physical_params.alpha, gmac, NULL); gfs_variables_swap (p, pmac); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_half_do, sim); gfs_centered_velocity_advection_diffusion (domain, FTT_DIMENSION, &sim->advection_params, gmac, sim->time.i > 0 || !gc ? gc : gmac, sim->physical_params.alpha); if (gc) { gfs_source_darcy_implicit (domain, sim->advection_params.dt); gfs_correct_centered_velocities (domain, FTT_DIMENSION, sim->time.i > 0 ? gc : gmac, -sim->advection_params.dt); /*inserted changes*/ gfs_post_projection (domain, por, FTT_DIMENSION); } else if (gfs_has_source_coriolis (domain)) { gfs_correct_centered_velocities (domain, FTT_DIMENSION, gmac, sim->advection_params.dt); gfs_source_darcy_implicit (domain, sim->advection_params.dt); gfs_correct_centered_velocities (domain, FTT_DIMENSION, gmac, -sim->advection_params.dt); /*inserted changes*/ gfs_post_projection (domain, por, FTT_DIMENSION); } gfs_domain_cell_traverse (domain, FTT_POST_ORDER, FTT_TRAVERSE_NON_LEAFS, -1, (FttCellTraverseFunc) gfs_cell_coarse_init, domain); gfs_simulation_adapt (sim); /*inserted changes */ gfs_approximate_projection_por (domain, por, &sim->approx_projection_params, sim->advection_params.dt, p, sim->physical_params.alpha, res, g, NULL); /*inserted changes */ sim->time.t = sim->tnext; sim->time.i++; gfs_simulation_set_timestep (sim); gfs_advance_tracers (sim, sim->advection_params.dt); gts_range_add_value (&domain->timestep, gfs_clock_elapsed (domain->timer) - tstart); gts_range_update (&domain->timestep); gts_range_add_value (&domain->size, gfs_domain_size (domain, FTT_TRAVERSE_LEAFS, -1)); gts_range_update (&domain->size); } gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_do, sim); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gts_object_destroy, NULL); for (c = 0; c < FTT_DIMENSION; c++) { gts_object_destroy (GTS_OBJECT (gmac[c])); if (sim->advection_params.gc) gts_object_destroy (GTS_OBJECT (g[c])); } }
extern field_t* field_read_gfs(const char* file) { int argc = 0; gfs_init(&argc,NULL); FILE *st = fopen(file,"r"); if (!st) { fprintf(stderr,"failed to open %s\n",file); return NULL; } GtsFile *fp = gts_file_new(st); GfsSimulation *sim = gfs_simulation_read(fp); if (!sim) { fprintf(stderr, "file %s not a valid simulation file\n" "line %d:%d: %s\n", file, fp->line, fp->pos, fp->error); return NULL; } gts_file_destroy(fp); fclose(st); GfsDomain *gdom = GFS_DOMAIN(sim); /* find the bounding box */ bbox_t *bb = NULL; gfs_domain_cell_traverse(gdom, FTT_PRE_ORDER, FTT_TRAVERSE_NON_LEAFS, 0, (FttCellTraverseFunc)ftt_bbox, &bb); if (!bb) { fprintf(stderr,"failed to determine bounding box\n"); return NULL; } #ifdef FRG_DEBUG fprintf(stdout, "%g %g %g %g\n", bb->x.min, bb->x.max, bb->y.min, bb->y.max); #endif /* tree depth and discretisation size */ int depth = gfs_domain_depth(gdom), nw = (int)(POW2(depth)*bbox_width(*bb)), nh = (int)(POW2(depth)*bbox_height(*bb)); /* shave bbox for node-aligned rather than pixel */ double shave = bbox_width(*bb)/(2.0*nw); bb->x.min += shave; bb->x.max -= shave; bb->y.min += shave; bb->y.max -= shave; /* create & intialise meshes */ bilinear_t* B[2]; int i; for (i=0 ; i<2 ; i++) { if ((B[i] = bilinear_new()) == NULL) return NULL; if (bilinear_dimension(nw,nh,*bb,B[i]) != ERROR_OK) return NULL; } ftts_t ftts; ftts.B = B; ftts.depth = depth; if ((ftts.u = gfs_variable_from_name(gdom->variables,"U")) == NULL) { fprintf(stderr,"no variable U\n"); return NULL; } if ((ftts.v = gfs_variable_from_name(gdom->variables,"V")) == NULL) { fprintf(stderr,"no variable V\n"); return NULL; } gfs_domain_cell_traverse(gdom, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc)ftt_sample, &ftts); /* clean up */ free(bb); gts_object_destroy(GTS_OBJECT(sim)); /* results */ field_t* F = malloc(sizeof(field_t)); if (!F) return NULL; F->u = B[0]; F->v = B[1]; F->k = NULL; return F; }
static void gfs_skew_symmetric_run (GfsSimulation * sim) { GfsVariable * p, * res = NULL, * gmac[FTT_DIMENSION]; GfsDomain * domain; GSList * i; domain = GFS_DOMAIN (sim); p = gfs_variable_from_name (domain->variables, "P"); g_assert (p); FttComponent c; for (c = 0; c < FTT_DIMENSION; c++) gmac[c] = gfs_temporary_variable (domain); gfs_variable_set_vector (gmac, FTT_DIMENSION); gfs_simulation_refine (sim); gfs_simulation_init (sim); i = domain->variables; while (i) { if (GFS_IS_VARIABLE_RESIDUAL (i->data)) res = i->data; i = i->next; } gfs_simulation_set_timestep (sim); GfsVariable ** u = gfs_domain_velocity (domain); GfsVariable ** velfaces = GFS_SKEW_SYMMETRIC(sim)->velfaces; GfsVariable ** velold = GFS_SKEW_SYMMETRIC(sim)->velold; FaceData fd = { velfaces, velold, u, p, &sim->advection_params.dt, GFS_SKEW_SYMMETRIC(sim)->beta}; if (sim->time.i == 0) { gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) reset_unold, &fd); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) get_face_values, &fd); gfs_mac_projection (domain, &sim->projection_params, sim->advection_params.dt/2., p, sim->physical_params.alpha, gmac, NULL); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) get_velfaces, &fd); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) initialize_unold, &fd); } while (sim->time.t < sim->time.end && sim->time.i < sim->time.iend) { gdouble tstart = gfs_clock_elapsed (domain->timer); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_do, sim); gfs_skew_symmetric_momentum (sim, &fd, gmac); gfs_mac_projection (domain, &sim->projection_params, sim->advection_params.dt/2., p, sim->physical_params.alpha, gmac, NULL); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) correct_face_velocity, NULL); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_half_do, sim); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) get_velfaces, &fd); gfs_domain_cell_traverse (domain, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttCellTraverseFunc) get_cell_values, &fd); gfs_domain_cell_traverse (domain, FTT_POST_ORDER, FTT_TRAVERSE_NON_LEAFS, -1, (FttCellTraverseFunc) gfs_cell_coarse_init, domain); gfs_simulation_adapt (sim); sim->time.t = sim->tnext; sim->time.i++; gfs_simulation_set_timestep (sim); gfs_advance_tracers (sim, sim->advection_params.dt); gts_range_add_value (&domain->timestep, gfs_clock_elapsed (domain->timer) - tstart); gts_range_update (&domain->timestep); gts_range_add_value (&domain->size, gfs_domain_size (domain, FTT_TRAVERSE_LEAFS, -1)); gts_range_update (&domain->size); } gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_do, sim); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gts_object_destroy, NULL); for (c = 0; c < FTT_DIMENSION; c++) gts_object_destroy (GTS_OBJECT (gmac[c])); }
static void wave_run (GfsSimulation * sim) { GfsDomain * domain = GFS_DOMAIN (sim); GfsWave * wave = GFS_WAVE (sim); SolidFluxParams par; par.div = gfs_variable_from_name (domain->variables, "P"); g_assert (par.div); par.p = &sim->advection_params; par.fv = gfs_temporary_variable (domain); gfs_simulation_refine (sim); gfs_simulation_init (sim); while (sim->time.t < sim->time.end && sim->time.i < sim->time.iend) { gdouble tstart = gfs_clock_elapsed (domain->timer); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_do, sim); /* get global timestep */ gfs_domain_face_traverse (domain, FTT_XYZ, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttFaceTraverseFunc) gfs_face_reset_normal_velocity, NULL); gfs_simulation_set_timestep (sim); gdouble dt = sim->advection_params.dt; gdouble g = sim->physical_params.g/sim->physical_params.L; gdouble tnext = sim->tnext; /* spatial advection */ guint ik, ith; for (ik = 0; ik < wave->nk; ik++) { FttVector cg; group_velocity (ik, 0, &cg, wave->ntheta, g); gfs_domain_face_traverse (domain, FTT_XYZ, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttFaceTraverseFunc) set_group_velocity, &cg); if (wave->alpha_s > 0.) { /* stability criterion for GSE diffusion */ gdouble cfl = sim->advection_params.cfl; sim->advection_params.cfl = MIN (cfl, 2./(4.*wave->alpha_s*M_PI/wave->ntheta)); /* fixme: this should be: sim->advection_params.cfl = MIN (cfl, sqrt(3.)/(wave->alpha_s*2.*M_PI/wave->ntheta)); */ gfs_simulation_set_timestep (sim); sim->advection_params.cfl = cfl; } else gfs_simulation_set_timestep (sim); /* subcycling */ guint n = rint (dt/sim->advection_params.dt); g_assert (fabs (sim->time.t + sim->advection_params.dt*n - tnext) < 1e-12); while (n--) { for (ith = 0; ith < wave->ntheta; ith++) { FttVector cg; group_velocity (ik, ith, &cg, wave->ntheta, g); gfs_domain_face_traverse (domain, FTT_XYZ, FTT_PRE_ORDER, FTT_TRAVERSE_LEAFS, -1, (FttFaceTraverseFunc) set_group_velocity, &cg); GfsVariable * t = GFS_WAVE (sim)->F[ik][ith]; sim->advection_params.v = t; gfs_domain_traverse_leaves (domain, (FttCellTraverseFunc) solid_flux, &par); gfs_tracer_advection_diffusion (domain, &sim->advection_params, NULL); sim->advection_params.fv = par.fv; gfs_domain_traverse_merged (domain, (GfsMergedTraverseFunc) gfs_advection_update, &sim->advection_params); if (wave->alpha_s > 0.) gse_alleviation_diffusion (domain, t, &cg, sim->advection_params.dt); gfs_domain_bc (domain, FTT_TRAVERSE_LEAFS, -1, t); gfs_domain_cell_traverse (domain, FTT_POST_ORDER, FTT_TRAVERSE_NON_LEAFS, -1, (FttCellTraverseFunc) t->fine_coarse, t); } gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) redo_some_events, sim); gfs_simulation_adapt (sim); } } sim->advection_params.dt = dt; /* source terms */ if (wave->source) (* wave->source) (wave); sim->time.t = sim->tnext = tnext; sim->time.i++; gts_range_add_value (&domain->timestep, gfs_clock_elapsed (domain->timer) - tstart); gts_range_update (&domain->timestep); gts_range_add_value (&domain->size, gfs_domain_size (domain, FTT_TRAVERSE_LEAFS, -1)); gts_range_update (&domain->size); } gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gfs_event_do, sim); gts_container_foreach (GTS_CONTAINER (sim->events), (GtsFunc) gts_object_destroy, NULL); gts_object_destroy (GTS_OBJECT (par.fv)); }
static void output_spectra_read (GtsObject ** o, GtsFile * fp) { (* GTS_OBJECT_CLASS (gfs_output_spectra_class ())->parent_class->read) (o, fp); if (fp->type == GTS_ERROR) return; GfsOutputSpectra * v = GFS_OUTPUT_SPECTRA (*o); if (fp->type != GTS_STRING) { gts_file_error (fp, "expecting a string (v)"); return; } GfsDomain * domain = GFS_DOMAIN (gfs_object_simulation (*o)); if (!(v->v = gfs_variable_from_name (domain->variables, fp->token->str))) { gts_file_error (fp, "unknown variable `%s'", fp->token->str); return; } gts_file_next_token (fp); GtsFileVariable var[] = { {GTS_DOUBLE, "x", TRUE, &v->pos.x}, {GTS_DOUBLE, "y", TRUE, &v->pos.y}, {GTS_DOUBLE, "z", TRUE, &v->pos.z}, {GTS_DOUBLE, "Lx", TRUE, &v->L.x}, {GTS_DOUBLE, "Ly", TRUE, &v->L.y}, {GTS_DOUBLE, "Lz", TRUE, &v->L.z}, {GTS_NONE} }; gts_file_assign_variables (fp, var); if (fp->type == GTS_ERROR) return; if (fp->type == GTS_INT) { v->level = atoi (fp->token->str); gts_file_next_token (fp); } else v->level = gfs_domain_depth (domain); guint i, j, k, size = 1; v->cgd = gfs_cartesian_grid_new (gfs_cartesian_grid_class ()); /* number of dims of the fft */ v->cgd->N = 3; v->Ndim = 0; k = 0; for (i = 0; i < 3; i++) { if ((&(v->L.x))[i] != 0) { v->Ndim++; v->dir[k] = i; k++; } } if (v->Ndim == 0) { gts_file_error (fp, "There must be at least one L component larger than 0"); return; } /* number of points in each direction */ v->cgd->n = g_malloc (3*sizeof (guint)); for (i = 0; i < 3; i++) { if ((&(v->L.x))[i] == 0 ) v->cgd->n[i] = 1; else v->cgd->n[i] = pow(2,v->level); size *= v->cgd->n[i]; } /* mesh coordinates */ v->cgd->x = g_malloc0 (3*sizeof (gdouble *)); for (i = 0; i < 3; i++) { v->cgd->x[i] = g_malloc (v->cgd->n[i]*sizeof (gdouble)); if (v->cgd->n[i] != 1) for (j = 0; j < v->cgd->n[i]; j++) v->cgd->x[i][j] = (&(v->pos.x))[i] + (&(v->L.x))[i]*(gdouble)j/((gdouble)(v->cgd->n[i]-1)) - 0.5; else v->cgd->x[i][0] = (&(v->pos.x))[i]; } /* memory data allocation */ v->cgd->v = g_malloc0( sizeof ( gdouble ) * 2*(size/2+1) ); }