123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401# 1 "src/owl/dense/owl_dense_matrix_generic.ml"(*
* OWL - OCaml Scientific and Engineering Computing
* Copyright (c) 2016-2020 Liang Wang <liang.wang@cl.cam.ac.uk>
*)openBigarrayopenOwl_ndarrayopenOwl_base_dense_commonincludeOwl_dense_ndarray_generic(* area function *)letshapex=letx_shape=shapexinassert(Array.lengthx_shape=2);x_shape.(0),x_shape.(1)(* creation functions *)letemptykmn=Owl_dense_ndarray_generic.emptyk[|m;n|]letcreatekmna=Owl_dense_ndarray_generic.createk[|m;n|]aletzeroskmn=Owl_dense_ndarray_generic.zerosk[|m;n|]letoneskmn=Owl_dense_ndarray_generic.onesk[|m;n|]letinitkmnf=Owl_dense_ndarray_generic.initk[|m;n|]fletinit_2dkmnf=letx=emptykmninlety=Bigarray.array2_of_genarrayxinfori=0tom-1doforj=0ton-1doBigarray.Array2.unsafe_setyij(fij)donedone;xleteyekn=Owl_dense_ndarray_generic.eyeknletsequentialk?a?stepmn=Owl_dense_ndarray_generic.sequentialk?a?step[|m;n|]letlinspacekabn=letx=Owl_dense_ndarray_generic.linspacekabninreshapex[|1;n|]letlogspacek?baseabn=letx=Owl_dense_ndarray_generic.logspacek?baseabninOwl_dense_ndarray_generic.reshapex[|1;n|]letdiagm?(k=0)v=letopenStdlibinletn=numelvinletu=reshapev[|n|]inletx=zeros(kindv)(n+absk)(n+absk)inleti,j=matchk<0with|true->absk,0|false->0,abskinfork=0ton-1dosetx[|i+k;j+k|](getu[|k|])done;xlettriu?(k=0)x=let_kind=kindxinletm,n=shapexinlety=zeros_kindmninletofs=refStdlib.(minn(max0k))inletlen=refStdlib.(max0(minn(n-k)))inletloops=Stdlib.(max0(minm(n-k)))infori=0toloops-1do_owl_copy_kind!len~ofsx:!ofs~incx:1~ofsy:!ofs~incy:1xy;ifi+k>=0then(ofs:=!ofs+n+1;len:=!len-1)elseofs:=!ofs+ndone;(* return the final upper triangular matrix *)ylettril?(k=0)x=let_kind=kindxinletm,n=shapexinlety=zeros_kindmninletrow_i=Stdlib.(minm(abs(min0k)))inletlen=refStdlib.(minn(max0k+1))inletofs=ref(row_i*n)infor_i=row_itom-1do_owl_copy_kind!len~ofsx:!ofs~incx:1~ofsy:!ofs~incy:1xy;ofs:=!ofs+n;if!len<nthenlen:=!len+1done;(* return the final lower triangular matrix *)yletsymmetric?(upper=true)x=let_kind=kindxinletm,n=shapexinOwl_exception.(check(m=n)(NOT_SQUARE[|m;n|]));lety=copyxinletofs=ref0inletincx,incy=matchupperwith|true->1,m|false->m,1infori=0tom-1do_owl_copy_kind(m-i)~ofsx:!ofs~incx~ofsy:!ofs~incyyy;ofs:=!ofs+n+1done;(* return the symmetric matrix *)yletbidiagonal?(upper=true)dvev=letm=numeldvinletn=numelevinassert(m-n=1);letx=zeros(kinddv)mminlet_x=Bigarray.array2_of_genarrayxinlet_dv=flattendv|>Bigarray.array1_of_genarrayinlet_ev=flattenev|>Bigarray.array1_of_genarrayinleti,j=matchupperwith|true->0,1|false->1,0infork=0tom-2do_x.{k,k}<-_dv.{k};_x.{k+i,k+j}<-_ev.{k}done;_x.{m-1,m-1}<-_dv.{m-1};xlethermitian?(upper=true)x=letm,n=shapexinOwl_exception.(check(m=n)(NOT_SQUARE[|m;n|]));lety=copyxinlet_y=flatteny|>Bigarray.array1_of_genarrayinlet_conj_op=_owl_conj(kindx)inletofs=ref0inletincx,incy=matchupperwith|true->1,m|false->m,1infori=0tom-1do(* copy and conjugate *)_conj_op(m-i)~ofsx:!ofs~incx~ofsy:!ofs~incyyy;(* set the imaginary part to zero by default. *)leta=_y.{!ofs}in_y.{!ofs}<-Complex.{re=a.re;im=0.};ofs:=!ofs+n+1done;(* return the symmetric matrix *)yletuniformk?a?bmn=Owl_dense_ndarray_generic.uniformk?a?b[|m;n|]letgaussiank?mu?sigmamn=Owl_dense_ndarray_generic.gaussiank?mu?sigma[|m;n|]letpoissonk~mumn=Owl_dense_ndarray_generic.poissonk~mu[|m;n|]letbernoullik?pmn=Owl_dense_ndarray_generic.bernoullik?p[|m;n|]letunit_basiskni=letx=Owl_dense_ndarray_generic.unit_basiskniinreshapex[|1;n|]lettoeplitz?cr=letc=matchcwith|Somec->c|None->conjrinletm=numelcinletn=numelrinOwl_dense_ndarray_generic.(setc[|0;0|](getr[|0;0|]));let_kind=kindrinletx=empty_kindmninletofs=ref0inletloops=Stdlib.minmninfori=0toloops-1do_owl_copy_kind(n-i)~ofsx:0~incx:1~ofsy:!ofs~incy:1rx;_owl_copy_kind(m-i)~ofsx:0~incx:1~ofsy:!ofs~incy:ncx;ofs:=!ofs+n+1done;xlethankel?rc=letm=numelcinletr=matchrwith|Somer->r|None->zeros(kindc)1minletn=numelrinOwl_dense_ndarray_generic.(setr[|0;0|](getc[|0;m-1|]));let_kind=kindrinletx=empty_kindmninletofs=ref((m-1)*n)inletloops=Stdlib.minmninfori=0toloops-1do_owl_copy_kind(n-i)~ofsx:0~incx:1~ofsy:!ofs~incy:1rx;_owl_copy_kind(m-i)~ofsx:i~incx:1~ofsy:i~incy:ncx;ofs:=!ofs-n+1done;xletkronab=letint_floorqr=int_of_float(float_of_intq/.float_of_intr)inletnra,nca=shapeainletnrb,ncb=shapebinletnrk,nck=nra*nrb,nca*ncbinletk=empty(kinda)nrknckinlet_mul_op=_mul_elt(kinda)inlet_a=Bigarray.array2_of_genarrayainlet_b=Bigarray.array2_of_genarraybinlet_k=Bigarray.array2_of_genarraykinfori=1tonrkdoforj=1tonckdoletai=int_floor(i-1)nrb+1inletaj=int_floor(j-1)ncb+1inletbi=((i-1)modnrb)+1inletbj=((j-1)modncb)+1in_k.{i-1,j-1}<-_mul_op_a.{ai-1,aj-1}_b.{bi-1,bj-1}donedone;k(* obtain properties *)letgetxij=Owl_dense_ndarray_generic.getx[|i;j|]letsetxija=Owl_dense_ndarray_generic.setx[|i;j|]aletswap_rowsxii'=letm,n=shapexinOwl_matrix._matrix_swap_rows(kindx)xmnii'letswap_colsxjj'=letm,n=shapexinOwl_matrix._matrix_swap_cols(kindx)xmnjj'lettransposex=letk=kindxinletm,n=shapexinlety=emptyknminOwl_matrix._matrix_transposekxy;yletctransposex=letk=kindxinletm,n=shapexinlety=emptyknminOwl_matrix._matrix_ctransposekxy;y(* iteration functions *)letiteri_rowsfx=fori=0torow_numx-1dofi(rowxi)doneletiter_rowsfx=iteri_rows(fun_y->fy)xletiter2i_rowsfxy=assert(row_numx=row_numy);iteri_rows(funiu->fiu(rowyi))xletiter2_rowsfxy=iter2i_rows(fun_uv->fuv)xylet_rowxi=(* get row i of x, but return as a column vector *)lety=slice_leftx[|i|]inreshapey[|col_numx;1|]letiteri_colsfx=lety=transposexinfori=0tocol_numx-1dofi(_rowyi)doneletiter_colsfx=iteri_cols(fun_y->fy)xletmapi_rowsfx=Array.init(row_numx)(funi->fi(rowxi))letmap_rowsfx=mapi_rows(fun_y->fy)xletmapi_colsfx=lety=transposexinArray.init(col_numx)(funi->fi(_rowyi))letmap_colsfx=mapi_cols(fun_y->fy)xletmapi_by_rowdfx=lety=empty(kindx)(row_numx)diniteri_rows(funiz->copy_row_to(fiz)yi)x;yletmap_by_rowdfx=mapi_by_rowd(fun_y->fy)xletmapi_by_coldfx=lety=empty(kindx)d(col_numx)initeri_cols(funjz->copy_col_to(fjz)yj)x;yletmap_by_coldfx=mapi_by_cold(fun_y->fy)xletmapi_at_rowfxi=letv=mapi(funjy->fjy)(rowxi)inlety=copyxincopy_row_tovyi;yletmap_at_rowfxi=mapi_at_row(fun_y->fy)xiletmapi_at_colfxj=letv=mapi(funiy->fiy)(colxj)inlety=copyxincopy_col_tovyj;yletmap_at_colfxj=mapi_at_col(fun_y->fy)xjletfilteri_rowsfx=lets=Owl_utils.Stack.make()initeri_rows(funiv->iffivthenOwl_utils.Stack.pushsi)x;Owl_utils.Stack.to_arraysletfilter_rowsfx=filteri_rows(fun_v->fv)xletfilteri_colsfx=lets=Owl_utils.Stack.make()initeri_cols(funiv->iffivthenOwl_utils.Stack.pushsi)x;Owl_utils.Stack.to_arraysletfilter_colsfx=filteri_cols(fun_v->fv)xlet_fold_basiciter_funfax=letr=refainiter_fun(funy->r:=f!ry)x;!rletfold_rowsfax=_fold_basiciter_rowsfaxletfold_colsfax=_fold_basiciter_colsfaxletiteri_2dfx=lety=array2_of_genarrayxinletm,n=shapexinfori=0tom-1doforj=0ton-1dofij(Array2.unsafe_getyij)donedoneletmapi_2dfx=lety=copyxinletz=array2_of_genarrayyinletm,n=shapeyinfori=0tom-1doforj=0ton-1doleta=Array2.unsafe_getzijinArray2.unsafe_setzij(fija)donedone;yletfilteri_2dfx=lets=Owl_utils.Stack.make()inlety=array2_of_genarrayxinletm,n=shapexinfori=0tom-1doforj=0ton-1doleta=Array2.unsafe_getyijiniffija=truethenOwl_utils.Stack.pushs(i,j)donedone;Owl_utils.Stack.to_arraysletfoldi_2d?axisfax=foldi_nd?axis(funkaccb->leti=k.(0)andj=k.(1)infijaccb)axletscani_2d?axisfx=scani_nd?axis(funkab->leti=k.(0)andj=k.(1)infijab)xletiter2i_2dfxy=assert(same_shapexy);letm,n=shapexinletx=array2_of_genarrayxinlety=array2_of_genarrayyinfori=0tom-1doforj=0ton-1doleta=Array2.unsafe_getxijinletb=Array2.unsafe_getyijinfijabdonedoneletmap2i_2dfxy=assert(same_shapexy);letm,n=shapexinletz=copyxinletx=array2_of_genarrayzinlety=array2_of_genarrayyinfori=0tom-1doforj=0ton-1doleta=Array2.unsafe_getxijinletb=Array2.unsafe_getyijinArray2.unsafe_setxij(fijab)donedone;zletsum_colsx=sum~axis:1xletsum_rowsx=sum~axis:0xletmean_colsx=mean~axis:1xletmean_rowsx=mean~axis:0xletmin_colsx=mapi_cols(funjv->letr,i=min_ivinr,i.(0),j)xletmin_rowsx=mapi_rows(funiv->letr,j=min_ivinr,i,j.(1))xletmax_colsx=mapi_cols(funjv->letr,i=max_ivinr,i.(0),j)xletmax_rowsx=mapi_rows(funiv->letr,j=max_ivinr,i,j.(1))xletmean'x=_mean_elt(kindx)(sum'x)(numelx)letadd_diagxa=letm,n=shapexinletm=Stdlib.minmninlety=copyxinlet_op=_add_elt(kindx)infori=0tom-1doletb=Owl_dense_ndarray_generic.getx[|i;i|]inOwl_dense_ndarray_generic.sety[|i;i|](_opba)done;y(* io operations *)letof_arraykxmn=letopenBigarrayinlety=Array1.of_arraykC_layoutx|>genarray_of_array1inOwl_dense_ndarray_generic.reshapey[|m;n|]letsave_txt?(sep="\t")?(append=false)~outx=letperm=0o666in(* will be AND'ed with user's umask *)letopen_flags=ifappendthen[Open_wronly;Open_creat;Open_append;Open_text]else[Open_wronly;Open_creat;Open_trunc;Open_text]inlet_op=Owl_utils.elt_to_str(kindx)inleth=open_out_genopen_flagspermoutiniter_rows(funy->iter(funz->Printf.fprintfh"%s%s"(_opz)sep)y;Printf.fprintfh"\n")x;close_outhletload_txt?(sep="\t")kf=let_op=Owl_utils.elt_of_strkinleth=open_infinlets=input_linehinletn=List.length(Str.split(Str.regexpsep)s)inletm=ref1in(* count lines in the input file *)(trywhiletruedoignore(input_lineh);m:=!m+1donewith|End_of_file->());letx=zerosk!mninseek_inh0;fori=0to!m-1dolets=Str.split(Str.regexpsep)(input_lineh)inList.iteri(funjy->setxij(_opy))sdone;close_inh;xletsemidefkn=letx=uniformknnindot(transposex)xletshuffle_rowsx=letm,_n=shapexinlety=copyxinfori=0tom-1doswap_rowsyi(Owl_stats.uniform_int_rvs~a:0~b:(m-1))done;yletshuffle_colsx=let_m,n=shapexinlety=copyxinfori=0ton-1doswap_colsyi(Owl_stats.uniform_int_rvs~a:0~b:(n-1))done;yletshufflex=x|>shuffle_rows|>shuffle_colsletmeshgridkxaxbyaybxnyn=letu=linspacekxaxbxninletv=linspacekyaybyninletx=map_by_rowxn(fun_->u)(emptykynxn)inlety=map_by_rowyn(fun_->v)(emptykxnyn)inx,transposeyletmeshupxy=letk=kindxinletxn=numelxinletyn=numelyinletx=map_by_rowxn(fun_->x)(emptykynxn)inlety=map_by_rowyn(fun_->y)(emptykxnyn)inx,transposey(* Hadamard Matrix *)let_hadamard_12=Array.mapfloat_of_int[|1;1;1;1;1;1;1;1;1;1;1;1;1;-1;1;-1;1;1;1;-1;-1;-1;1;-1;1;-1;-1;1;-1;1;1;1;-1;-1;-1;1;1;1;-1;-1;1;-1;1;1;1;-1;-1;-1;1;-1;1;-1;-1;1;-1;1;1;1;-1;-1;1;-1;-1;1;-1;-1;1;-1;1;1;1;-1;1;-1;-1;-1;1;-1;-1;1;-1;1;1;1;1;1;-1;-1;-1;1;-1;-1;1;-1;1;1;1;1;1;-1;-1;-1;1;-1;-1;1;-1;1;1;1;1;1;-1;-1;-1;1;-1;-1;1;-1;1;-1;1;1;1;-1;-1;-1;1;-1;-1;1;1;1;-1;1;1;1;-1;-1;-1;1;-1;-1|]let_hadamard_20=Array.mapfloat_of_int[|1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;1;-1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;1;1;-1;-1;1;1;-1;-1;-1;-1;1;-1;1;-1;1;1;1;1;-1;-1|]lethadamardkn=(* function to build up hadamard matrix recursively *)letrec_make_hadamard(cp_op:('a,'b)Owl_core_types.owl_arr_op18)(neg_op:('a,'b)Owl_core_types.owl_arr_op18)lennbasex=iflen=basethen()else(letlen'=len/2in_make_hadamardcp_opneg_oplen'nbasex;letofsx=ref0infor_i=0tolen'-1doletx1_ofs=!ofsx+len'inletx2_ofs=!ofsx+(len'*n)inletx3_ofs=x2_ofs+len'incp_oplen'~ofsx:!ofsx~incx:1~ofsy:x1_ofs~incy:1xx;cp_oplen'~ofsx:!ofsx~incx:1~ofsy:x2_ofs~incy:1xx;cp_oplen'~ofsx:!ofsx~incx:1~ofsy:x3_ofs~incy:1xx;(* negate the bottom right block *)neg_oplen'~ofsx:x3_ofs~incx:1~ofsy:x3_ofs~incy:1xx;ofsx:=!ofsx+ndone)in(* function to convert the pre-calculated hadamard array into type k *)let_float_array_to_k:typeab.(a,b)kind->floatarray->aarray=funka->matchkwith|Bigarray.Float32->a|Bigarray.Float64->a|Bigarray.Complex32->Array.map(funb->Complex.{re=b;im=0.})a|Bigarray.Complex64->Array.map(funb->Complex.{re=b;im=0.})a|_->failwith"Owl_dense_matrix_generic.hadamard"in(* start building, only deal with pow2 of n, n/12, n/20. *)letx=emptyknninletcp_op=_owl_copykinletneg_op=_owl_negkinifOwl_maths.is_pow2nthen(Owl_dense_ndarray_generic.setx[|0;0|](Owl_const.onek);_make_hadamardcp_opneg_opnn1x;x)elseifOwl_maths.is_pow2(n/12)&&Stdlib.(nmod12)=0then(lety=_float_array_to_kk_hadamard_12inlety=of_arrayky1212inlet_area=area001111incopy_area_toy_areax_area;_make_hadamardcp_opneg_opnn12x;x)elseifOwl_maths.is_pow2(n/20)&&Stdlib.(nmod20)=0then(lety=_float_array_to_kk_hadamard_20inlety=of_arrayky2020inlet_area=area001919incopy_area_toy_areax_area;_make_hadamardcp_opneg_opnn20x;x)elsefailwith"Owl_dense_matrix_generic:hadamard"letmagickn=leterror()=lets=Printf.sprintf"magic requires n >= 3, whereas input n = %i"ninOwl_exception.INVALID_ARGUMENTsinOwl_exception.verify(n>=3)error;leta0=Owl_const.zerokinleta1=Owl_const.onekinlet_magic_oddna=letx=zerosknninleti=ref0inletj=ref(n/2)inletac=ref(float_of_inta|>_float_typ_eltk)inletm=(n*n)+a-1|>float_of_int|>_float_typ_eltkinwhile!ac<=mdoifgetx!i!j=a0then(setx!i!j!ac;ac:=_add_eltk!aca1)else(leti'=(!i-1+n)modninletj'=(!j+1)modninifgetxi'j'=a0then(i:=i';j:=j')elsei:=(!i+1)modn)done;xinlet_magic_doubly_evenn=letx=zerosknninlet_seq_incxi0j0i1j1=fori=i0toi1doletac=ref((n*i)+j0+1|>float_of_int|>_float_typ_eltk)inforj=j0toj1dosetxij!ac;ac:=_add_eltk!aca1donedoneinlet_seq_decxi0j0i1j1=letac=ref(n*n|>float_of_int|>_float_typ_eltk)infori=i0toi1doforj=j0toj1doifgetxij=a0thensetxij!ac;ac:=_sub_eltk!aca1donedoneinletm=n/4in_seq_incx00(m-1)(m-1);_seq_incx0(3*m)(m-1)((4*m)-1);_seq_incxmm((3*m)-1)((3*m)-1);_seq_incx(3*m)0((4*m)-1)(m-1);_seq_incx(3*m)(3*m)((4*m)-1)((4*m)-1);_seq_decx00(n-1)(n-1);xinlet_magic_singly_evenn=letm=n/2inletm2=m*minletxa=_magic_oddm1inletxb=_magic_oddm(m2+1)inletxc=_magic_oddm((2*m2)+1)inletxd=_magic_oddm((3*m2)+1)inletm3=m/2inletxa'=concat_horizontal(get_slice[[];[0;m3-1]]xd)(get_slice[[];[m3;-1]]xa)inletxd'=concat_horizontal(get_slice[[];[0;m3-1]]xa)(get_slice[[];[m3;-1]]xd)inletxb'=ifm3-1=0thenxbelseconcat_horizontal(get_slice[[];[0;m-m3]]xb)(get_slice[[];[1-m3;-1]]xc)inletxc'=ifm3-1=0thenxcelseconcat_horizontal(get_slice[[];[0;m-m3]]xc)(get_slice[[];[1-m3;-1]]xb)insetxa'm30(getxam30);setxa'm3m3(getxdm3m3);setxd'm30(getxdm30);setxd'm3m3(getxam3m3);concat_horizontal(concat_verticalxa'xd')(concat_verticalxc'xb')in(* n is odd *)ifOwl_maths.is_oddnthen_magic_oddn1(* n is doubly even *)elseifnmod4=0then_magic_doubly_evenn(* n is singly even *)else_magic_singly_evennletmax_pool?paddingxkernelstride=letm,n=shapexinletx=Owl_dense_ndarray_generic.reshapex[|1;m;n;1|]inlety=Owl_dense_ndarray_generic.max_pool2d?paddingxkernelstrideinlets=Owl_dense_ndarray_generic.shapeyinletm,n=s.(1),s.(2)inOwl_dense_ndarray_generic.reshapey[|m;n|]letavg_pool?paddingxkernelstride=letm,n=shapexinletx=Owl_dense_ndarray_generic.reshapex[|1;m;n;1|]inlety=Owl_dense_ndarray_generic.avg_pool2d?paddingxkernelstrideinlets=Owl_dense_ndarray_generic.shapeyinletm,n=s.(1),s.(2)inOwl_dense_ndarray_generic.reshapey[|m;n|]letcov?b~a=leta=matchbwith|Someb->letna=numelainletnb=numelbinassert(na=nb);leta=reshapea[|na;1|]inletb=reshapeb[|nb;1|]inconcat_horizontalab|None->ainletmu=mean_rowsainleta=subamuinleta'=ctransposeainletc=dota'ainletn=row_numa-1|>Stdlib.max1|>float_of_int|>_float_typ_elt(kinda)indiv_scalarcnletmat2gray?amin?amaxx=letamin=matchaminwith|Somea->a|None->min'xinletamax=matchamaxwith|Somea->a|None->max'xinletx=clip_by_value~amin~amaxxinletx=sub_scalarxaminindiv_scalarx(_sub_elt(kindx)amaxamin)(* end here *)