OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
lapacke_cbbcsd.c
Go to the documentation of this file.
1/*****************************************************************************
2 Copyright (c) 2014, Intel Corp.
3 All rights reserved.
4
5 Redistribution and use in source and binary forms, with or without
6 modification, are permitted provided that the following conditions are met:
7
8 * Redistributions of source code must retain the above copyright notice,
9 this list of conditions and the following disclaimer.
10 * Redistributions in binary form must reproduce the above copyright
11 notice, this list of conditions and the following disclaimer in the
12 documentation and/or other materials provided with the distribution.
13 * Neither the name of Intel Corporation nor the names of its contributors
14 may be used to endorse or promote products derived from this software
15 without specific prior written permission.
16
17 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF
27 THE POSSIBILITY OF SUCH DAMAGE.
28*****************************************************************************
29* Contents: Native high-level C interface to LAPACK function cbbcsd
30* Author: Intel Corporation
31*****************************************************************************/
32
33#include "lapacke_utils.h"
34
35lapack_int LAPACKE_cbbcsd( int matrix_layout, char jobu1, char jobu2,
36 char jobv1t, char jobv2t, char trans, lapack_int m,
37 lapack_int p, lapack_int q, float* theta, float* phi,
42 float* b11d, float* b11e, float* b12d, float* b12e,
43 float* b21d, float* b21e, float* b22d, float* b22e )
44{
45 lapack_int info = 0;
46 lapack_int lrwork = -1;
47 float* rwork = NULL;
48 float rwork_query;
49 int lapack_layout;
50 if( matrix_layout != LAPACK_COL_MAJOR && matrix_layout != LAPACK_ROW_MAJOR ) {
51 LAPACKE_xerbla( "LAPACKE_cbbcsd", -1 );
52 return -1;
53 }
54 if( LAPACKE_lsame( trans, 'n' ) && matrix_layout == LAPACK_COL_MAJOR ) {
55 lapack_layout = LAPACK_COL_MAJOR;
56 } else {
57 lapack_layout = LAPACK_ROW_MAJOR;
58 }
59#ifndef LAPACK_DISABLE_NAN_CHECK
60 if( LAPACKE_get_nancheck() ) {
61 /* Optionally check input matrices for NaNs */
62 if( LAPACKE_s_nancheck( q-1, phi, 1 ) ) {
63 return -11;
64 }
65 if( LAPACKE_s_nancheck( q, theta, 1 ) ) {
66 return -10;
67 }
68 if( LAPACKE_lsame( jobu1, 'y' ) ) {
69 if( LAPACKE_cge_nancheck( lapack_layout, p, p, u1, ldu1 ) ) {
70 return -12;
71 }
72 }
73 if( LAPACKE_lsame( jobu2, 'y' ) ) {
74 if( LAPACKE_cge_nancheck( lapack_layout, m-p, m-p, u2, ldu2 ) ) {
75 return -14;
76 }
77 }
78 if( LAPACKE_lsame( jobv1t, 'y' ) ) {
79 if( LAPACKE_cge_nancheck( lapack_layout, q, q, v1t, ldv1t ) ) {
80 return -16;
81 }
82 }
83 if( LAPACKE_lsame( jobv2t, 'y' ) ) {
84 if( LAPACKE_cge_nancheck( lapack_layout, m-q, m-q, v2t, ldv2t ) ) {
85 return -18;
86 }
87 }
88 }
89#endif
90 /* Query optimal working array(s) size */
91 info = LAPACKE_cbbcsd_work( matrix_layout, jobu1, jobu2, jobv1t, jobv2t,
92 trans, m, p, q, theta, phi, u1, ldu1, u2, ldu2,
93 v1t, ldv1t, v2t, ldv2t, b11d, b11e, b12d, b12e,
94 b21d, b21e, b22d, b22e, &rwork_query, lrwork );
95 if( info != 0 ) {
96 goto exit_level_0;
97 }
98 lrwork = (lapack_int)rwork_query;
99 /* Allocate memory for work arrays */
100 rwork = (float*)LAPACKE_malloc( sizeof(float) * lrwork );
101 if( rwork == NULL ) {
103 goto exit_level_0;
104 }
105 /* Call middle-level interface */
106 info = LAPACKE_cbbcsd_work( matrix_layout, jobu1, jobu2, jobv1t, jobv2t,
107 trans, m, p, q, theta, phi, u1, ldu1, u2, ldu2,
108 v1t, ldv1t, v2t, ldv2t, b11d, b11e, b12d, b12e,
109 b21d, b21e, b22d, b22e, rwork, lrwork );
110 /* Release memory and exit */
111 LAPACKE_free( rwork );
112exit_level_0:
113 if( info == LAPACK_WORK_MEMORY_ERROR ) {
114 LAPACKE_xerbla( "LAPACKE_cbbcsd", info );
115 }
116 return info;
117}
#define lapack_int
Definition lapack.h:83
#define lapack_complex_float
Definition lapack.h:45
#define LAPACK_WORK_MEMORY_ERROR
Definition lapacke.h:55
#define LAPACK_COL_MAJOR
Definition lapacke.h:53
#define LAPACKE_free(p)
Definition lapacke.h:46
#define LAPACK_ROW_MAJOR
Definition lapacke.h:52
lapack_int LAPACKE_cbbcsd_work(int matrix_layout, char jobu1, char jobu2, char jobv1t, char jobv2t, char trans, lapack_int m, lapack_int p, lapack_int q, float *theta, float *phi, lapack_complex_float *u1, lapack_int ldu1, lapack_complex_float *u2, lapack_int ldu2, lapack_complex_float *v1t, lapack_int ldv1t, lapack_complex_float *v2t, lapack_int ldv2t, float *b11d, float *b11e, float *b12d, float *b12e, float *b21d, float *b21e, float *b22d, float *b22e, float *rwork, lapack_int lrwork)
int LAPACKE_get_nancheck(void)
#define LAPACKE_malloc(size)
Definition lapacke.h:43
lapack_int LAPACKE_cbbcsd(int matrix_layout, char jobu1, char jobu2, char jobv1t, char jobv2t, char trans, lapack_int m, lapack_int p, lapack_int q, float *theta, float *phi, lapack_complex_float *u1, lapack_int ldu1, lapack_complex_float *u2, lapack_int ldu2, lapack_complex_float *v1t, lapack_int ldv1t, lapack_complex_float *v2t, lapack_int ldv2t, float *b11d, float *b11e, float *b12d, float *b12e, float *b21d, float *b21e, float *b22d, float *b22e)
lapack_logical LAPACKE_lsame(char ca, char cb)
void LAPACKE_xerbla(const char *name, lapack_int info)
lapack_logical LAPACKE_s_nancheck(lapack_int n, const float *x, lapack_int incx)
lapack_logical LAPACKE_cge_nancheck(int matrix_layout, lapack_int m, lapack_int n, const lapack_complex_float *a, lapack_int lda)