static void p4est_coarsen_both (p4est_t * p4est, int coarsen_recursive, p4est_coarsen_t coarsen_fn, p4est_init_t init_fn) { int success; p4est_t *copy; copy = p4est_copy (p4est, 1); p4est_coarsen_old (copy, coarsen_recursive, coarsen_fn, init_fn); coarsen_callback_count = 0; p4est_coarsen_ext (p4est, coarsen_recursive, 1, coarsen_fn, init_fn, NULL); SC_CHECK_ABORT (coarsen_recursive || coarsen_callback_count == (int) p4est->local_num_quadrants, "Coarsen count"); success = p4est_is_equal (p4est, copy, 1); SC_CHECK_ABORT (success, "Coarsen mismatch"); p4est_destroy (copy); }
int main (int argc, char **argv) { int mpirank, mpisize; int mpiret; sc_MPI_Comm mpicomm; p4est_t *p4est; p4est_connectivity_t *connectivity; mpiret = sc_MPI_Init (&argc, &argv); SC_CHECK_MPI (mpiret); mpicomm = sc_MPI_COMM_WORLD; mpiret = sc_MPI_Comm_size (mpicomm, &mpisize); SC_CHECK_MPI (mpiret); mpiret = sc_MPI_Comm_rank (mpicomm, &mpirank); SC_CHECK_MPI (mpiret); sc_init (mpicomm, 1, 1, NULL, SC_LP_DEFAULT); p4est_init (NULL, SC_LP_DEFAULT); /* create connectivity and forest structures */ #ifdef P4_TO_P8 connectivity = p8est_connectivity_new_rotcubes (); #else connectivity = p4est_connectivity_new_star (); #endif p4est = p4est_new_ext (mpicomm, connectivity, 15, 0, 0, 1, NULL, NULL); p4est_refine_ext (p4est, 1, P4EST_QMAXLEVEL, refine_fn, NULL, replace_fn); p4est_coarsen_ext (p4est, 1, 0, coarsen_fn, NULL, replace_fn); p4est_balance_ext (p4est, P4EST_CONNECT_FULL, NULL, replace_fn); p4est_destroy (p4est); p4est_connectivity_destroy (connectivity); sc_finalize (); mpiret = sc_MPI_Finalize (); SC_CHECK_MPI (mpiret); return 0; }
int p4est_wrap_adapt (p4est_wrap_t * pp) { int changed; #ifdef P4EST_ENABLE_DEBUG p4est_locidx_t jl, local_num; #endif p4est_gloidx_t global_num; p4est_t *p4est = pp->p4est; P4EST_ASSERT (!pp->hollow); P4EST_ASSERT (pp->coarsen_delay >= 0); P4EST_ASSERT (pp->mesh != NULL); P4EST_ASSERT (pp->ghost != NULL); P4EST_ASSERT (pp->mesh_aux == NULL); P4EST_ASSERT (pp->ghost_aux == NULL); P4EST_ASSERT (pp->match_aux == 0); P4EST_ASSERT (pp->temp_flags == NULL); P4EST_ASSERT (pp->num_refine_flags >= 0 && pp->num_refine_flags <= p4est->local_num_quadrants); /* This allocation is optimistic when not all refine requests are honored */ pp->temp_flags = P4EST_ALLOC_ZERO (uint8_t, p4est->local_num_quadrants + (P4EST_CHILDREN - 1) * pp->num_refine_flags); /* Execute refinement */ pp->inside_counter = pp->num_replaced = 0; #ifdef P4EST_ENABLE_DEBUG local_num = p4est->local_num_quadrants; #endif global_num = p4est->global_num_quadrants; p4est_refine_ext (p4est, 0, -1, refine_callback, NULL, replace_on_refine); P4EST_ASSERT (pp->inside_counter == local_num); P4EST_ASSERT (p4est->local_num_quadrants - local_num == pp->num_replaced * (P4EST_CHILDREN - 1)); changed = global_num != p4est->global_num_quadrants; /* Execute coarsening */ pp->inside_counter = pp->num_replaced = 0; #ifdef P4EST_ENABLE_DEBUG local_num = p4est->local_num_quadrants; #endif global_num = p4est->global_num_quadrants; p4est_coarsen_ext (p4est, 0, 1, coarsen_callback, NULL, pp->coarsen_delay ? replace_on_coarsen : pp->replace_fn); P4EST_ASSERT (pp->inside_counter == local_num); P4EST_ASSERT (local_num - p4est->local_num_quadrants == pp->num_replaced * (P4EST_CHILDREN - 1)); changed = changed || global_num != p4est->global_num_quadrants; /* Free temporary flags */ P4EST_FREE (pp->temp_flags); pp->temp_flags = NULL; /* Only if refinement and/or coarsening happened do we need to balance */ if (changed) { P4EST_FREE (pp->flags); p4est_balance_ext (p4est, pp->btype, NULL, pp->coarsen_delay ? replace_on_balance : pp->replace_fn); pp->flags = P4EST_ALLOC_ZERO (uint8_t, p4est->local_num_quadrants); pp->ghost_aux = p4est_ghost_new (p4est, pp->btype); pp->mesh_aux = p4est_mesh_new_ext (p4est, pp->ghost_aux, 1, 1, pp->btype); pp->match_aux = 1; } #ifdef P4EST_ENABLE_DEBUG else { for (jl = 0; jl < p4est->local_num_quadrants; ++jl) { P4EST_ASSERT (pp->flags[jl] == 0); } } #endif pp->num_refine_flags = 0; return changed; }
/** Timestep the advection problem. * * Update the state, refine, repartition, and write the solution to file. * * \param [in,out] p4est the forest, whose state is updated * \param [in] time the end time */ static void step3_timestep (p4est_t * p4est, double time) { double t = 0.; double dt = 0.; int i; step3_data_t *ghost_data; step3_ctx_t *ctx = (step3_ctx_t *) p4est->user_pointer; int refine_period = ctx->refine_period; int repartition_period = ctx->repartition_period; int write_period = ctx->write_period; int recursive = 0; int allowed_level = P4EST_QMAXLEVEL; int allowcoarsening = 1; int callbackorphans = 0; int mpiret; double orig_max_err = ctx->max_err; double umax, global_umax; p4est_ghost_t *ghost; /* create the ghost quadrants */ ghost = p4est_ghost_new (p4est, P4EST_CONNECT_FULL); /* create space for storing the ghost data */ ghost_data = P4EST_ALLOC (step3_data_t, ghost->ghosts.elem_count); /* synchronize the ghost data */ p4est_ghost_exchange_data (p4est, ghost, ghost_data); /* initialize du/dx estimates */ p4est_iterate (p4est, ghost, (void *) ghost_data, /* pass in ghost data that we just exchanged */ step3_reset_derivatives, /* blank the previously calculated derivatives */ step3_minmod_estimate, /* compute the minmod estimate of each cell's derivative */ #ifdef P4_TO_P8 NULL, /* there is no callback for the edges between quadrants */ #endif NULL); /* there is no callback for the corners between quadrants */ for (t = 0., i = 0; t < time; t += dt, i++) { P4EST_GLOBAL_PRODUCTIONF ("time %f\n", t); /* refine */ if (!(i % refine_period)) { if (i) { /* compute umax */ umax = 0.; /* initialize derivative estimates */ p4est_iterate (p4est, NULL, (void *) &umax, /* pass in ghost data that we just exchanged */ step3_compute_max, /* blank the previously calculated derivatives */ NULL, /* there is no callback for the faces between quadrants */ #ifdef P4_TO_P8 NULL, /* there is no callback for the edges between quadrants */ #endif NULL); /* there is no callback for the corners between quadrants */ mpiret = sc_MPI_Allreduce (&umax, &global_umax, 1, sc_MPI_DOUBLE, sc_MPI_MAX, p4est->mpicomm); SC_CHECK_MPI (mpiret); ctx->max_err = orig_max_err * global_umax; P4EST_GLOBAL_PRODUCTIONF ("u_max %f\n", global_umax); /* adapt */ p4est_refine_ext (p4est, recursive, allowed_level, step3_refine_err_estimate, NULL, step3_replace_quads); p4est_coarsen_ext (p4est, recursive, callbackorphans, step3_coarsen_err_estimate, NULL, step3_replace_quads); p4est_balance_ext (p4est, P4EST_CONNECT_FACE, NULL, step3_replace_quads); p4est_ghost_destroy (ghost); P4EST_FREE (ghost_data); ghost = NULL; ghost_data = NULL; } dt = step3_get_timestep (p4est); } /* repartition */ if (i && !(i % repartition_period)) { p4est_partition (p4est, allowcoarsening, NULL); if (ghost) { p4est_ghost_destroy (ghost); P4EST_FREE (ghost_data); ghost = NULL; ghost_data = NULL; } } /* write out solution */ if (!(i % write_period)) { step3_write_solution (p4est, i); } /* synchronize the ghost data */ if (!ghost) { ghost = p4est_ghost_new (p4est, P4EST_CONNECT_FULL); ghost_data = P4EST_ALLOC (step3_data_t, ghost->ghosts.elem_count); p4est_ghost_exchange_data (p4est, ghost, ghost_data); } /* compute du/dt */ /* *INDENT-OFF* */ p4est_iterate (p4est, /* the forest */ ghost, /* the ghost layer */ (void *) ghost_data, /* the synchronized ghost data */ step3_quad_divergence, /* callback to compute each quad's interior contribution to du/dt */ step3_upwind_flux, /* callback to compute each quads' faces' contributions to du/du */ #ifdef P4_TO_P8 NULL, /* there is no callback for the edges between quadrants */ #endif NULL); /* there is no callback for the corners between quadrants */ /* *INDENT-ON* */ /* update u */ p4est_iterate (p4est, NULL, /* ghosts are not needed for this loop */ (void *) &dt, /* pass in dt */ step3_timestep_update, /* update each cell */ NULL, /* there is no callback for the faces between quadrants */ #ifdef P4_TO_P8 NULL, /* there is no callback for the edges between quadrants */ #endif NULL); /* there is no callback for the corners between quadrants */ /* synchronize the ghost data */ p4est_ghost_exchange_data (p4est, ghost, ghost_data); /* update du/dx estimate */ p4est_iterate (p4est, ghost, (void *) ghost_data, /* pass in ghost data that we just exchanged */ step3_reset_derivatives, /* blank the previously calculated derivatives */ step3_minmod_estimate, /* compute the minmod estimate of each cell's derivative */ #ifdef P4_TO_P8 NULL, /* there is no callback for the edges between quadrants */ #endif NULL); /* there is no callback for the corners between quadrants */ } P4EST_FREE (ghost_data); p4est_ghost_destroy (ghost); }