int main (int argc, char **argv) { #ifdef VSG_HAVE_MPI VsgPRTreeParallelConfig pconfig = {{NULL,}}; #endif VsgVector2d lbound = {-1., -1.}; VsgVector2d ubound = {1., 1.}; VsgPRTree2d *prtree; AranSolver2d *solver; int ret = 0; guint i; GTimer *timer = NULL; #ifdef VSG_HAVE_MPI MPI_Init (&argc, &argv); MPI_Comm_size (MPI_COMM_WORLD, &sz); MPI_Comm_rank (MPI_COMM_WORLD, &rk); #endif aran_init(); parse_args (argc, argv); #ifdef VSG_HAVE_MPI pconfig.communicator = MPI_COMM_WORLD; pconfig.point = point_accum_vtable; aran_development2d_vtable_init (&pconfig.node_data, 0, order); #endif points = g_ptr_array_new (); if (check) check_points = g_malloc0 (np * sizeof (PointAccum)); prtree = vsg_prtree2d_new_full (&lbound, &ubound, (VsgPoint2dLocFunc) vsg_vector2d_vector2d_locfunc, (VsgPoint2dDistFunc) vsg_vector2d_dist, (VsgRegion2dLocFunc) NULL, maxbox); aran_binomial_require (2*order); solver = aran_solver2d_new (prtree, ARAN_TYPE_DEVELOPMENT2D, aran_development2d_new (0, order), (AranZeroFunc) aran_development2d_set_zero); #ifdef VSG_HAVE_MPI aran_solver2d_set_parallel (solver, &pconfig); #endif if (virtual_maxbox != 0) aran_solver2d_set_nf_isleaf (solver, _nf_isleaf_virtual_maxbox, &virtual_maxbox); aran_solver2d_set_functions (solver, (AranParticle2ParticleFunc2d) p2p, (AranParticle2MultipoleFunc2d) p2m, (AranMultipole2MultipoleFunc2d) aran_development2d_m2m, (AranMultipole2LocalFunc2d) aran_development2d_m2l, (AranLocal2LocalFunc2d) aran_development2d_l2l, (AranLocal2ParticleFunc2d)l2p); if (semifar_threshold < G_MAXUINT) { aran_solver2d_set_functions_full (solver, (AranParticle2ParticleFunc2d) p2p, (AranParticle2MultipoleFunc2d) p2m, (AranMultipole2MultipoleFunc2d) aran_development2d_m2m, (AranMultipole2LocalFunc2d) aran_development2d_m2l, (AranLocal2LocalFunc2d) aran_development2d_l2l, (AranLocal2ParticleFunc2d) l2p, (AranParticle2LocalFunc2d) p2l, (AranMultipole2ParticleFunc2d) m2p, semifar_threshold); if (semifar_threshold == 0) { PointAccum p1 = {{0.1, 0.1}, 0.1, 0., 0}; PointAccum p2 = {{-0.1, -0.1}, 0.1, 0., 1}; /* compute operators timings to be able to compute optimal solver parameters */ aran_solver2d_profile_operators (solver, (AranParticleInitFunc2d) point_accum_clear_accum, &p1, &p2); /* alternatively, we could get timings from profile databases */ /* aran_profile_db_read_file ("./profiledb-newtonfield3.ini", NULL); */ /* aran_solver2d_db_profile_operators (solver, (gdouble) order); */ } } if (_hilbert) { /* configure for hilbert curve order traversal */ aran_solver2d_set_children_order_hilbert (solver); } _distribution (points, solver); if (_verbose) { g_printerr ("%d : solve begin\n", rk); #ifdef VSG_HAVE_MPI MPI_Barrier (MPI_COMM_WORLD); #endif timer = g_timer_new (); } aran_solver2d_solve (solver); if (_verbose) { #ifdef VSG_HAVE_MPI MPI_Barrier (MPI_COMM_WORLD); #endif g_printerr ("%d : solve ok elapsed=%f seconds\n", rk, g_timer_elapsed (timer, NULL)); g_timer_destroy (timer); } if (check) { if (sz == 1) { for (i=0; i<np; i++) { PointAccum *pi = &check_points[i]; guint j; for (j=0; j<np; j++) { if (i != j) { PointAccum *pj = &check_points[j]; gcomplex128 zd_m_zs = (pi->vector.x + I*pi->vector.y) - (pj->vector.x + I*pj->vector.y); pi->accum += 1./zd_m_zs * pj->density; } } } } else check_parallel_points (solver); aran_solver2d_foreach_point (solver, (GFunc) check_point_accum, &ret); if (_verbose) g_printerr ("%d : max err = %e\n", rk, maxerr); g_free (check_points); } aran_solver2d_free (solver); #ifdef VSG_HAVE_MPI aran_development2d_vtable_clear (&pconfig.node_data); #endif /* destroy the points */ g_ptr_array_foreach (points, empty_array, NULL); g_ptr_array_free (points, TRUE); #ifdef VSG_HAVE_MPI MPI_Finalize (); #endif return ret; }
int main (int argc, char **argv) { VsgVector2d lbound = {-1., -1.}; VsgVector2d ubound = {1., 1.}; PointAccum *points[N] = {0}; VsgPRTree2d *prtree; AranSolver2d *solver; int ret = 0; guint i; aran_init(); parse_args (argc, argv); prtree = vsg_prtree2d_new_full (&lbound, &ubound, (VsgPoint2dLocFunc) vsg_vector2d_vector2d_locfunc, (VsgPoint2dDistFunc) vsg_vector2d_dist, (VsgRegion2dLocFunc) NULL, maxbox); aran_binomial_require (2*order); solver = aran_solver2d_new (prtree, ARAN_TYPE_DEVELOPMENT2D, aran_development2d_new (0, order), (AranZeroFunc) aran_development2d_set_zero); aran_solver2d_set_functions (solver, (AranParticle2ParticleFunc2d) p2p, (AranParticle2MultipoleFunc2d) p2m, (AranMultipole2MultipoleFunc2d) aran_development2d_m2m, (AranMultipole2LocalFunc2d) aran_development2d_m2l, (AranLocal2LocalFunc2d) aran_development2d_l2l, (AranLocal2ParticleFunc2d)l2p); _distribution (points, solver); aran_solver2d_solve (solver); /* vsg_prtree2d_write (prtree, stderr); */ aran_solver2d_free (solver); if (check) { for (i=0; i<np; i++) { guint j; gcomplex128 sum = 0.; gcomplex128 err; for (j=0; j<np; j++) { if (i != j) { gcomplex128 zd_m_zs = (points[i]->vector.x + I*points[i]->vector.y) - (points[j]->vector.x + I*points[j]->vector.y); sum += 1./zd_m_zs * points[j]->density; } } err = (points[i]->accum - sum) / MAX(cabs (points[i]->accum), cabs (sum)); if (cabs (err) > err_lim) { g_printerr ("Error: pt%u (%e,%e) -> ", i, creal (err), cimag (err)); vsg_vector2d_write (&points[i]->vector, stderr); g_printerr ("\n"); ret ++; } } } for (i=0; i<np; i++) { g_free (points[i]); } return ret; }
void check_parallel_points (AranSolver2d *tocheck) { VsgVector2d lbound = {-1., -1.}; VsgVector2d ubound = {1., 1.}; VsgPRTree2d *prtree; AranSolver2d *solver; int i; prtree = vsg_prtree2d_new_full (&lbound, &ubound, (VsgPoint2dLocFunc) vsg_vector2d_vector2d_locfunc, (VsgPoint2dDistFunc) vsg_vector2d_dist, (VsgRegion2dLocFunc) NULL, maxbox); aran_binomial_require (2*order); solver = aran_solver2d_new (prtree, ARAN_TYPE_DEVELOPMENT2D, aran_development2d_new (0, order), (AranZeroFunc) aran_development2d_set_zero); aran_solver2d_set_functions (solver, (AranParticle2ParticleFunc2d) p2p, (AranParticle2MultipoleFunc2d) p2m, (AranMultipole2MultipoleFunc2d) aran_development2d_m2m, (AranMultipole2LocalFunc2d) aran_development2d_m2l, (AranLocal2LocalFunc2d) aran_development2d_l2l, (AranLocal2ParticleFunc2d)l2p); if (semifar_threshold < G_MAXUINT) { /* if optimal threshold was requested, we need to compare with the same value */ if (semifar_threshold == 0) aran_solver2d_get_functions_full (tocheck, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &semifar_threshold); aran_solver2d_set_functions_full (solver, (AranParticle2ParticleFunc2d) p2p, (AranParticle2MultipoleFunc2d) p2m, (AranMultipole2MultipoleFunc2d) aran_development2d_m2m, (AranMultipole2LocalFunc2d) aran_development2d_m2l, (AranLocal2LocalFunc2d) aran_development2d_l2l, (AranLocal2ParticleFunc2d) l2p, (AranParticle2LocalFunc2d) p2l, (AranMultipole2ParticleFunc2d) m2p, semifar_threshold); } if (_hilbert) { /* configure for hilbert curve order traversal */ aran_solver2d_set_children_order_hilbert (solver); } for (i=0; i<np; i++) { aran_solver2d_insert_point (solver, &check_points[i]); } aran_solver2d_solve (solver); aran_solver2d_free (solver); }