Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

FEAT: Added additional tools for Fortran calling codes #1154

Merged
merged 3 commits into from
Jan 19, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
26 changes: 26 additions & 0 deletions include/cantera/base/NoExitLogger.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,26 @@
//! @file NoExitLogger.h

// This file is part of Cantera. See License.txt in the top-level directory or
// at https://cantera.org/license.txt for license and copyright information.

#ifndef NOEXITLOGGER_H
#define NOEXITLOGGER_H

#include <string>
#include "cantera/base/logger.h"

namespace Cantera {
/// Logger that doesn't exit when an error is thrown.
/// @ingroup textlogs
class NoExitLogger : public Logger {
public:
NoExitLogger() {}
virtual ~NoExitLogger() {}

virtual void error(const std::string& msg)
{
std::cerr << msg << std::endl;
}
};
}
#endif
13 changes: 13 additions & 0 deletions samples/f90/demo.f90
Original file line number Diff line number Diff line change
Expand Up @@ -58,6 +58,7 @@ subroutine demo(gas, MAXSP, MAXRXNS)

double precision q(MAXRXNS), qf(MAXRXNS), qr(MAXRXNS)
double precision diff(MAXSP)
double precision molar_cp(MAXSP), molar_h(MAXSP), g_rt(MAXSP)

character*80 eq
character*20 name
Expand Down Expand Up @@ -119,6 +120,18 @@ subroutine demo(gas, MAXSP, MAXRXNS)
40 format(' ',a20,e14.5,' m2/s')
end do

! Thermodynamic properties
CALL getPartialMolarCp(gas, molar_cp)
CALL getPartialMolarEnthalpies(gas, molar_h)
CALL getGibbs_RT(gas, g_rt)

write(*,*) 'Species molar cp, molar enthalpy, normalized Gibbs free energy'
do k = 1, nsp
call getSpeciesName(gas, k, name)
write(*,50) name, molar_cp(k), molar_h(k), g_rt(k)
50 format(' ',a20,g14.5,' J/kmol-K',g14.5,' J/kmol',g14.5,'-')
end do

return

end subroutine demo
16 changes: 16 additions & 0 deletions src/fortran/cantera.f90
Original file line number Diff line number Diff line change
Expand Up @@ -111,6 +111,10 @@ MODULE CANTERA
MODULE PROCEDURE ctfunc_getCanteraError
END INTERFACE getCanteraError

INTERFACE turnOffExitOnError
MODULE PROCEDURE ctfunc_TurnOffExitOnError
END INTERFACE turnOffExitOnError

INTERFACE getCp_R
MODULE PROCEDURE ctthermo_getCp_R
END INTERFACE getCp_R
Expand All @@ -131,6 +135,10 @@ MODULE CANTERA
MODULE PROCEDURE ctthermo_getEnthalpies_RT
END INTERFACE getEnthalpies_RT

INTERFACE getGibbs_RT
MODULE PROCEDURE ctthermo_getGibbs_RT
END INTERFACE getGibbs_RT

INTERFACE getEntropies_R
MODULE PROCEDURE ctthermo_getEntropies_R
END INTERFACE getEntropies_R
Expand Down Expand Up @@ -183,6 +191,14 @@ MODULE CANTERA
MODULE PROCEDURE ctthermo_getPartialMolarIntEnerg_R
END INTERFACE getPartialMolarIntEnergies

INTERFACE getPartialMolarCp
MODULE PROCEDURE ctthermo_getPartialMolarCp
END INTERFACE getPartialMolarCp

INTERFACE getPartialMolarEnthalpies
MODULE PROCEDURE ctthermo_getPartialMolarEnthalpies
END INTERFACE getPartialMolarEnthalpies

INTERFACE getReactionString
MODULE PROCEDURE ctkin_getReactionString
END INTERFACE getReactionString
Expand Down
5 changes: 5 additions & 0 deletions src/fortran/cantera_funcs.f90
Original file line number Diff line number Diff line change
Expand Up @@ -86,6 +86,11 @@ subroutine ctfunc_getCanteraError(buf)
ierr = ctgetCanteraError(buf)
end subroutine ctfunc_getCanteraError

subroutine ctfunc_turnOffExitOnError()
implicit none
call ctturnOffExitOnError()
end subroutine ctfunc_turnOffExitOnError
tpg2114 marked this conversation as resolved.
Show resolved Hide resolved

subroutine ctfunc_addCanteraDirectory(self, buf)
implicit none
type(phase_t), intent(inout) :: self
Expand Down
20 changes: 20 additions & 0 deletions src/fortran/cantera_thermo.f90
Original file line number Diff line number Diff line change
Expand Up @@ -545,6 +545,13 @@ subroutine ctthermo_getEnthalpies_RT(self, h_rt)
self%err = th_getenthalpies_rt(self%thermo_id, h_rt)
end subroutine ctthermo_getenthalpies_rt

subroutine ctthermo_getGibbs_RT(self, grt)
implicit none
type(phase_t), intent(inout) :: self
double precision, intent(out) :: grt(self%nsp)
self%err = th_getgibbs_rt(self%thermo_id, grt)
end subroutine ctthermo_getGibbs_RT

subroutine ctthermo_getEntropies_R(self, s_r)
implicit none
type(phase_t), intent(inout) :: self
Expand All @@ -567,4 +574,17 @@ subroutine ctthermo_getPartialMolarIntEnerg_R(self, ie)
self%err = th_getpartialmolarintenergies_r(self%thermo_id, ie)
end subroutine ctthermo_getPartialMolarIntEnerg_R

subroutine ctthermo_getPartialMolarCp(self, cpbar)
implicit none
type(phase_t), intent(inout) :: self
double precision, intent(out) :: cpbar(self%nsp)
self%err = th_getpartialmolarcp(self%thermo_id, cpbar)
end subroutine ctthermo_getPartialMolarCp

subroutine ctthermo_getPartialMolarEnthalpies(self, hbar)
implicit none
type(phase_t), intent(inout) :: self
double precision, intent(out) :: hbar(self%nsp)
self%err = th_getpartialmolarenthalpies(self%thermo_id, hbar)
end subroutine ctthermo_getPartialMolarEnthalpies
end module cantera_thermo
Loading