Path: blob/main/api/t8_fortran_interface/t8_fortran_interface.h
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/** \file t8_fortran_interface.h23* In this file we provide a basic Fortran interface24* for some functions of t8code.25* Mostly, the C functions here are wrappers around more complex26* t8code function.27* We only export a minimum of the actual t8code functionality28* to Fortran.29*/3031#ifndef T8_FORTRAN_INTERFACE_H32#define T8_FORTRAN_INTERFACE_H3334#include <t8.h>35#include <t8_cmesh.h>36#include <t8_forest/t8_forest_general.h>37#include <t8_forest/t8_forest_geometrical.h>3839typedef int (*t8_fortran_adapt_coordinate_callback) (double x, double y, double z, int is_family);4041/* A fallback type if t8code is not built with MPI */42typedef43#if T8_ENABLE_MPI44MPI_Fint45#else46int47#endif48MPI_T8_Fint;4950T8_EXTERN_C_BEGIN ();5152/** Initialize sc and t8code with SC_MPI_COMM_WORLD communicator53* and SC_LP_DEFAULT logging.54* This call is equivalent to55* sc_init (comm, 1, 1, NULL, SC_LP_ESSENTIAL);56* t8_init (SC_LP_DEFAULT);57* \param [in] comm The MPI communicator to use.58*/59void60t8_fortran_init_all (sc_MPI_Comm *comm);6162/** Finalize sc. This wraps sc_finalize in order to have consistent63* naming with t8_fortran_init_all.64*/65void66t8_fortran_finalize ();6768/** Commit cmesh. This wraps cmesh_commit in order to use the dereferenced communicator.69* \param [in, out] cmesh Cmesh to commit70* \param [in] Ccomm Pointer to a C MPI communicator.71*/72void73t8_fortran_cmesh_commit (t8_cmesh_t cmesh, sc_MPI_Comm *comm);7475/** This function calls t8_cmesh_set_join_by_stash with connectivity = NULL.76* \param[in,out] cmesh Pointer to a t8code cmesh object. If set to NULL this argument is ignored.77* \param[in] do_both_directions Compute the connectivity from both neighboring sides.78* Takes much longer to compute.79*80* \warning This routine might be too expensive for very large meshes. In this case,81* consider to use a fully featured mesh generator.82*83* \note This routine does not detect periodic boundaries.84*/85void86t8_fortran_cmesh_set_join_by_stash_noConn (t8_cmesh_t cmesh, const int do_both_directions);8788/** Translate a fortran MPI communicator into a C MPI communicator89* and return a pointer to it.90* \param [in] Fcomm Fortran MPI Communicator91* \return Pointer to the corresponding C MPI communicator.92* \note This function allocated memory for the new C MPI communicator.93* Call \ref t8_fortran_MPI_Comm_delete to free this memory.94* \note t8code needs to be configured with MPI support to be able to use95* this function.96*/97sc_MPI_Comm *98t8_fortran_MPI_Comm_new (MPI_T8_Fint Fcomm);99100/** Free the memory of a C MPI Communicator pointer that was created101* with \ref t8_fortran_MPI_Comm_new.102* \param [in] Ccomm Pointer to a C MPI communicator.103*/104void105t8_fortran_MPI_Comm_delete (sc_MPI_Comm *Ccomm);106107/** Wraps t8_cmesh_new_periodic_tri, passing the MPI communicator as pointer instead of by value108* \param [in] Ccomm Pointer to a C MPI communicator.109* \param [out] Example cmesh110*/111t8_cmesh_t112t8_cmesh_new_periodic_tri_wrap (sc_MPI_Comm *Ccomm);113114/** Wraps \ref t8_forest_new_uniform with the default scheme as scheme115* and passes MPI communicator as pointer instead of by value.116* Build a uniformly refined forest on a coarse mesh.117* \param [in] cmesh A coarse mesh.118* \param [in] level An initial uniform refinement level.119* \param [in] do_face_ghost If true, a layer of ghost elements is created for the forest.120* \param [in] comm MPI communicator to use.121* \return A uniform forest with coarse mesh \a cmesh, eclass_scheme122* \a scheme and refinement level \a level.123*/124t8_forest_t125t8_forest_new_uniform_default (t8_cmesh_t cmesh, int level, int do_face_ghost, sc_MPI_Comm *comm);126127/**128* \param [in, out] forest The forest129* \param [in] recursive A flag specifying whether adaptation is to be done recursively130* or not. If the value is zero, adaptation is not recursive131* and it is recursive otherwise.132* \param [in] callback A pointer to a user defined function. t8code will never touch the function.133*/134t8_forest_t135t8_forest_adapt_by_coordinates (t8_forest_t forest, int recursive, t8_fortran_adapt_coordinate_callback callback);136137/** Log a message on the root rank with priority SC_LP_PRODUCTION.138* \param [in] string String to log.139*/140void141t8_global_productionf_noargs (const char *string);142143T8_EXTERN_C_END ();144145#endif /* !T8_FORTRAN_INTERFACE_H */146147148