Path: blob/main/api/t8_fortran_interface/t8_fortran_interface.c
504 views
/*1This file is part of t8code.2t8code is a C library to manage a collection (a forest) of multiple3connected adaptive space-trees of general element classes in parallel.45Copyright (C) 2025 the developers67t8code is free software; you can redistribute it and/or modify8it under the terms of the GNU General Public License as published by9the Free Software Foundation; either version 2 of the License, or10(at your option) any later version.1112t8code is distributed in the hope that it will be useful,13but WITHOUT ANY WARRANTY; without even the implied warranty of14MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the15GNU General Public License for more details.1617You should have received a copy of the GNU General Public License18along with t8code; if not, write to the Free Software Foundation, Inc.,1951 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.20*/2122#include <t8_fortran_interface.h>23#include <t8_forest/t8_forest_general.h>24#include <t8_forest/t8_forest_geometrical.h>25#include <t8_cmesh/t8_cmesh_examples.h>26#include <t8_cmesh/t8_cmesh_helpers.h>27#include <t8_schemes/t8_scheme.h>28#include <t8_schemes/t8_default/t8_default_c_interface.h>2930void31t8_fortran_init_all_ (sc_MPI_Comm *comm)32{33T8_ASSERT (comm != NULL);34/* Initialize sc */35sc_init (*comm, 1, 1, NULL, SC_LP_DEFAULT);36/* Initialize t8code */37t8_init (SC_LP_DEFAULT);38}3940/* Wrapper around sc_finalize */41void42t8_fortran_finalize ()43{44sc_finalize ();45}4647void48t8_fortran_cmesh_commit (t8_cmesh_t cmesh, sc_MPI_Comm *comm)49{50t8_cmesh_commit (cmesh, *comm);51}5253void54t8_fortran_cmesh_set_join_by_stash_noConn (t8_cmesh_t cmesh, const int do_both_directions)55{56t8_cmesh_set_join_by_stash (cmesh, NULL, do_both_directions);57}5859void60t8_fortran_init_all (sc_MPI_Comm *comm)61{62int rank;6364t8_fortran_init_all_ (comm);65if (*comm != sc_MPI_COMM_NULL) {66sc_MPI_Comm_rank (*comm, &rank);67t8_debugf ("rank = %i\n", rank);68}69}7071void72t8_fortran_init_all_noMPI ()73{74sc_MPI_Comm commnull = sc_MPI_COMM_NULL;75t8_fortran_init_all (&commnull);76}7778/* Build C MPI comm from Fortran MPI Comm. */79sc_MPI_Comm *80t8_fortran_MPI_Comm_new (MPI_T8_Fint Fcomm)81{82#if !T8_ENABLE_MPI83SC_ABORT ("t8code was not configured with MPI support.");84return NULL;85#endif86/* We use malloc instead of T8_ALLOC since t8code may not be initialized87* yet. */88sc_MPI_Comm *Ccomm = (sc_MPI_Comm *) malloc (sizeof (*Ccomm));89#if T8_ENABLE_MPI90/* If configured with MPI, transform the Fortran communicator handle to a C handle */91*Ccomm = MPI_Comm_f2c (Fcomm);92#else93/* In case it is not configured with MPI, set the communicator to NULL as a fallback */94*Ccomm = sc_MPI_COMM_NULL;95#endif96t8_debugf ("Created comm %lu\n", (long unsigned) Ccomm);97return Ccomm;98}99100/* Delete C MPI Comm. */101void102t8_fortran_MPI_Comm_delete (sc_MPI_Comm *Ccomm)103{104#if !T8_ENABLE_MPI105SC_ABORT ("t8code was not configured with MPI support.");106#endif107free (Ccomm);108}109110t8_cmesh_t111t8_cmesh_new_periodic_tri_wrap (sc_MPI_Comm *Ccomm)112{113return t8_cmesh_new_periodic_tri (*Ccomm);114}115116t8_forest_t117t8_forest_new_uniform_default (t8_cmesh_t cmesh, int level, int do_face_ghost, sc_MPI_Comm *comm)118{119const t8_scheme_c *default_scheme = t8_scheme_new_default ();120121T8_ASSERT (comm != NULL);122return t8_forest_new_uniform (cmesh, default_scheme, level, do_face_ghost, *comm);123}124125int126t8_fortran_adapt_by_coordinates_callback (t8_forest_t forest, t8_forest_t forest_from, t8_locidx_t which_tree,127const t8_eclass_t tree_class,128__attribute__ ((unused)) t8_locidx_t lelement_id, const t8_scheme_c *scheme,129const int is_family, const int num_elements, t8_element_t *elements[])130{131t8_fortran_adapt_coordinate_callback callback132= (t8_fortran_adapt_coordinate_callback) t8_forest_get_user_function (forest);133double midpoint[3];134t8_forest_element_centroid (forest_from, which_tree, elements[0], midpoint);135t8_debugf ("Coord: %.2f\n", midpoint[0]);136int ret = callback (midpoint[0], midpoint[1], midpoint[2], num_elements > 0);137138/* Coarsen if a family was given and return value is negative. */139if (is_family) {140/* The elements form a family */141T8_ASSERT (t8_elements_are_family (scheme, tree_class, elements));142/* Build the parent. */143t8_element_t *parent;144t8_element_new (scheme, tree_class, 1, &parent);145t8_element_get_parent (scheme, tree_class, elements[0], parent);146/* Get the coordinates of the parent. */147t8_forest_element_centroid (forest_from, which_tree, parent, midpoint);148149ret = callback (midpoint[0], midpoint[1], midpoint[2], 1);150}151else {152/* The elements do not form a family. */153/* Get the coordinates of the first element and call callback */154t8_forest_element_centroid (forest_from, which_tree, elements[0], midpoint);155ret = callback (midpoint[0], midpoint[1], midpoint[2], 0);156T8_ASSERT (ret >= 0);157}158return ret;159}160161t8_forest_t162t8_forest_adapt_by_coordinates (t8_forest_t forest, int recursive, t8_fortran_adapt_coordinate_callback callback)163{164t8_forest_t forest_new;165166T8_ASSERT (t8_forest_is_committed (forest));167T8_ASSERT (callback != NULL);168169/* Initialize new forest */170t8_forest_init (&forest_new);171/* Set the callback as user data */172t8_forest_set_user_function (forest_new, (t8_generic_function_pointer) callback);173/* Call set adapt */174t8_forest_set_adapt (forest_new, forest, t8_fortran_adapt_by_coordinates_callback, recursive);175/* Commit the forest */176t8_forest_commit (forest_new);177return forest_new;178}179180void181t8_global_productionf_noargs (const char *string)182{183t8_global_productionf ("%s", string);184}185186187