1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175
|
/*
Copyright (C) 2021 Fredrik Johansson
This file is part of Calcium.
Calcium is free software: you can redistribute it and/or modify it under
the terms of the GNU Lesser General Public License (LGPL) as published
by the Free Software Foundation; either version 2.1 of the License, or
(at your option) any later version. See <http://www.gnu.org/licenses/>.
*/
#include "ca.h"
void
ca_asin_special(ca_t res, const ca_t x, ca_ctx_t ctx)
{
if (ca_check_is_signed_inf(x, ctx) == T_TRUE)
{
ca_t s;
ca_init(s, ctx);
/* -i*csgn(i*z)*inf */
ca_i(s, ctx);
ca_mul(res, x, s, ctx);
ca_csgn(res, res, ctx);
ca_mul(res, res, s, ctx);
ca_neg(res, res, ctx);
ca_pos_inf(s, ctx);
ca_mul(res, res, s, ctx);
ca_clear(s, ctx);
return;
}
if (ca_check_is_uinf(x, ctx) == T_TRUE || ca_check_is_undefined(x, ctx) == T_TRUE)
{
ca_set(res, x, ctx);
return;
}
ca_unknown(res, ctx);
return;
}
static int
_ca_asin_rational(ca_t res, const ca_t x, ca_ctx_t ctx)
{
qqbar_t v;
slong p;
ulong q;
int success;
qqbar_init(v);
/* todo: rule out non-sines more quickly */
if (ca_get_qqbar(v, x, ctx) && qqbar_asin_pi(&p, &q, v))
{
ca_pi(res, ctx);
ca_mul_si(res, res, p, ctx);
ca_div_ui(res, res, q, ctx);
success = 1;
}
else
{
success = 0;
}
qqbar_clear(v);
return success;
}
void
ca_asin_logarithm(ca_t res, const ca_t x, ca_ctx_t ctx)
{
ca_t t, u;
if (CA_IS_SPECIAL(x))
{
ca_asin_special(res, x, ctx);
return;
}
if (_ca_asin_rational(res, x, ctx))
return;
/* asin(x) = -i log(ix + sqrt(1-x^2)) */
ca_init(t, ctx);
ca_init(u, ctx);
ca_sqr(t, x, ctx);
ca_ui_sub(t, 1, t, ctx);
ca_sqrt(t, t, ctx);
ca_i(u, ctx);
ca_mul(u, u, x, ctx);
ca_add(t, t, u, ctx);
ca_log(t, t, ctx);
ca_i(u, ctx);
ca_mul(res, t, u, ctx);
ca_neg(res, res, ctx);
ca_clear(t, ctx);
ca_clear(u, ctx);
}
void
ca_asin_direct(ca_t res, const ca_t x, ca_ctx_t ctx)
{
if (CA_IS_SPECIAL(x))
{
ca_asin_special(res, x, ctx);
return;
}
if (_ca_asin_rational(res, x, ctx))
return;
/* todo: csgn normalization, reflection...? */
_ca_function_fx(res, CA_Asin, x, ctx);
}
void
ca_asin(ca_t res, const ca_t x, ca_ctx_t ctx)
{
if (ctx->options[CA_OPT_TRIG_FORM] == CA_TRIG_EXPONENTIAL)
{
ca_asin_logarithm(res, x, ctx);
}
/* todo:
else if (ctx->options[CA_OPT_TRIG_FORM] == CA_TRIG_TANGENT)
{
ca_asin_arctangent(res, x, ctx);
}
*/
else
{
ca_asin_direct(res, x, ctx);
}
}
void
ca_acos_logarithm(ca_t res, const ca_t x, ca_ctx_t ctx)
{
ca_t t;
ca_init(t, ctx);
ca_pi(t, ctx);
ca_div_ui(t, t, 2, ctx);
ca_asin_logarithm(res, x, ctx);
ca_sub(res, t, res, ctx);
ca_clear(t, ctx);
}
void
ca_acos_direct(ca_t res, const ca_t x, ca_ctx_t ctx)
{
ca_t t;
ca_init(t, ctx);
ca_pi(t, ctx);
ca_div_ui(t, t, 2, ctx);
ca_asin_direct(res, x, ctx);
ca_sub(res, t, res, ctx);
ca_clear(t, ctx);
}
void
ca_acos(ca_t res, const ca_t x, ca_ctx_t ctx)
{
ca_t t;
ca_init(t, ctx);
ca_pi(t, ctx);
ca_div_ui(t, t, 2, ctx);
ca_asin(res, x, ctx);
ca_sub(res, t, res, ctx);
ca_clear(t, ctx);
}
|