BioFSharp


Accessing online databases

BioFSharp contains a set of readers that facilitate the access to different biological online resources. This documentation aims to give an introduction for them.

UniProt's Proteins REST API online access

📂View module documentation
The Proteins REST API provides access to key biological data from UniProt and data from Large Scale Studies data mapped to UniProt. The services provide sequence feature annotations from UniProtKB, variation data from UniProtKB and mapped from Large Scale data sources (1 Genomes, ExAC and COSMIC), proteomics data mapped from Large Scale sources (PeptideAtlas, MaxQB and EPD) and genome coordinate mappings.

In this tutorial I want to show you how to access data from UniProt and give examples for what might be done with it. As a first step we will retreive data about the Copper-transporting ATPase 2 of different species.

1: 
2: 
3: 
4: 
let mouse = EbiAPI.UniProteinDB.getProteinSeqFeature "Q64446" |> Async.RunSynchronously
let rat = EbiAPI.UniProteinDB.getProteinSeqFeature "Q64535"   |> Async.RunSynchronously
let human = EbiAPI.UniProteinDB.getProteinSeqFeature "P35670" |> Async.RunSynchronously
let sheep = EbiAPI.UniProteinDB.getProteinSeqFeature "Q9XT50" |> Async.RunSynchronously
No value has been returned

What we want to do now is getting an measure of phylogenetic relationship. We will do this by aligning the peptide-sequences with each other. For this we prepare the data and the alignment function.

1: 
2: 
3: 
let proteinInfos = [|human;mouse;rat;sheep|]
let sequences = proteinInfos |> Array.map ((fun x -> x.Sequence) >> BioArray.ofAminoAcidSymbolString)
let names = proteinInfos |> Array.map (fun x -> x.EntryName)

 1: 
 2: 
 3: 
 4: 
 5: 
 6: 
 7: 
 8: 
 9: 
10: 
11: 
12: 
13: 
14: 
15: 
16: 
17: 
18: 
19: 
20: 
21: 
///Creates the pairwise global alignment of two amino acid sequences based on the BLOSUM62 scoring matrix
let align a b =
    let scoring = ScoringMatrix.getScoringMatrixAminoAcid ScoringMatrix.BLOSUM62
    let costs : PairwiseAlignment.Costs<AminoAcidSymbols.AminoAcidSymbol> = 
        {Open = -4;
        Continuation = -3;
        Similarity = scoring}
    PairwiseAlignment.NeedlemanWunsch.runAminoAcidSymbol costs a b

///Aligns two amino acid sequences and writes their relative alignment score into a square matrix
let createAlignmentMatrix (sequences : BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> []) = 
        let mat = Array2D.zeroCreate 4 4
        for i = 0 to 3 do
            for j = i to 3 do 
                let alignment = align sequences.[i] sequences.[j]
                mat.[i,j] <- float alignment.MetaData
                mat.[j,i] <- float alignment.MetaData
        let mat' = Array2D.map (fun x ->  x / (Array2D.maxBy id mat)) mat
        for i = 0 to 3 do
            mat'.[i,i] <- 1.
        mat'

Let's see what we get..
1: 
2: 
3: 
4: 
let heatMap =
    createAlignmentMatrix sequences
    |> Array2D.toJaggedArray
    |> fun x -> Chart.Heatmap(data = x,RowNames = names, ColNames = names,Colorscale = StyleParam.Colorscale.Custom([0.,"#ffffff";1.,"#44546A"]),Name = "Relative pairwise alignment score")
No value has been returned

We can see that the Copper-transporting ATPase 2 of mouse and rat are the most similar to each other, while the one of the sheep is relatively different.

But not only the protein names and their sequences can be stored in UniProt. For these proteins we also get information about the cellular location of each of their amino acids. Why not see if the hydrophobic amino acids are located in the membrane after all?
First we design a type to model the location more nicely. Afterwards we create a function for mapping the protein information to just the location information of all amino acids in one protein.

 1: 
 2: 
 3: 
 4: 
 5: 
 6: 
 7: 
 8: 
 9: 
10: 
11: 
12: 
13: 
14: 
15: 
16: 
17: 
18: 
19: 
20: 
21: 
type Location = 
    |Intracellular
    |Transmembrane
    |Extracellular
///Maps protein information to information of location of amino acids of this protein
let getLocationInfo (proteinInfo:EbiAPI.ProteinsAPIschema.ProteinFeatureInfo) = 
    proteinInfo.Features
    |> Array.choose (fun x -> 
        match x.Type with
        | "TOPO_DOM" | "TRANSMEM" -> Some (x.Begin,x.End,x.Description)
        | _ -> None)
    |> Array.collect (fun (s,e,info) -> 
                        Array.init 
                            ((int e) - (int s)+1)
                            (fun i -> 
                                match info with
                                | _ when info = "Cytoplasmic"   -> Intracellular
                                | _ when info = "Helical"       -> Transmembrane
                                | _ when info = "Extracellular" -> Extracellular
                            )
                        )

We also create ourselves a function for retreiving hydrophobicity of all amino acids of our protein and an additional helper.

 1: 
 2: 
 3: 
 4: 
 5: 
 6: 
 7: 
 8: 
 9: 
10: 
11: 
12: 
///Evaluates the hydrophobicity of all amino acids of a sequence
let getHydrophobicity a =
    let hydro = AminoProperties.initGetAminoProperty AminoProperties.HydrophobicityIndex
    AminoProperties.ofWindowedBioArray 10 hydro a

///Creates the squares depicting the cell membrane
let createBackgroundShape (i,(location : Location)) =
    let i = float i 
    match location with 
    | Transmembrane -> 
        Some (Shape.init (ShapeType = StyleParam.ShapeType.Rectangle, X0 = (i-0.5), X1 = (i+0.5),Y0 = 0.,Y1  = 2.,Opacity = 0.05,Fillcolor = Colors.toWebColor Colors.Table.Office.orange))
    | _ -> None
We pick the mouse protein for this analysis.
1: 
2: 
3: 
let hydrophobicityArr = 
    getHydrophobicity sequences.[1]
let mouseLocInfo = mouse |> getLocationInfo

Now let's see what we get. The colored areas depict the cell membrane. The line graph represents the amino acids of the protein with their given hydrophobicity.

1: 
2: 
3: 
4: 
5: 
let shapes = 
    Chart.Line ((Array.indexed hydrophobicityArr),Color = Colors.toHex false Colors.Table.Office.blue)
    |> Chart.withShapes(mouseLocInfo |> Seq.indexed |> Seq.choose createBackgroundShape)
    |> Chart.withX_AxisStyle("Sequenceindex")
    |> Chart.withY_AxisStyle("Hydrophobicity")
No value has been returned

As you can see the hydrophobic areas overlap with the membrane areas.

namespace BioFSharp
namespace BioFSharp.BioDB
Multiple items
namespace FSharp

--------------------
namespace Microsoft.FSharp
namespace FSharp.Plotly
namespace FSharpAux
namespace BioFSharp.Algorithm
module Colors

from FSharpAux
module Formula

from BioFSharp
module Table

from BioFSharp.Formula
val mouse : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo
module EbiAPI

from BioFSharp.BioDB
Multiple items
type UniProteinDB =
  new : unit -> UniProteinDB
  static member getAntigen : accession:string -> Async<ProteinFeatureInfo>
  static member getGenomicCoordinates : accession:string -> Async<GnEntry>
  static member getProteinSeqFeature : accession:string * ?categories:string [] * ?types:string [] -> Async<ProteinFeatureInfo>
  static member getProteomicsPeptide : accession:string -> Async<ProteinFeatureInfo>
  static member getUniProtEntry : accession:string -> Async<Entry>
  static member searchAntigens : ?offset:int * ?size:int * ?accession:string * ?antigenSequence:string * ?antigenID:string * ?ensemblIDs:string * ?matchscore:int -> Async<ProteinFeatureInfo []>
  static member searchEntries : ?offset:int * ?size:int * ?accession:string * ?reviewed:string * ?isoforms:int * ?goterms:string * ?keywords:string * ?ec:string * ?gene:string * ?protein:string * ?organism:string * ?taxid:string * ?pubmed:string -> Async<Entry []>
  static member searchGenomicCoordinates : ?offset:int * ?size:int * ?accession:string * ?chromosome:string * ?ensembl:string * ?gene:string * ?protein:string * ?taxid:string * ?location:string -> Async<GnEntry []>
  static member searchNaturalVariants : ?offset:int * ?size:int * ?accession:string * ?taxid:string * ?sourcetype:string * ?consequencetype:string * ?wildtype:string * ?alternativesequence:string * ?location:string * ?disease:string * ?omim:string * ?evidence:string -> Async<ProteinFeatureInfo []>
  ...

--------------------
new : unit -> EbiAPI.UniProteinDB
static member EbiAPI.UniProteinDB.getProteinSeqFeature : accession:string * ?categories:string [] * ?types:string [] -> Async<EbiAPI.ProteinsAPIschema.ProteinFeatureInfo>
Multiple items
type Async =
  static member AsBeginEnd : computation:('Arg -> Async<'T>) -> ('Arg * AsyncCallback * obj -> IAsyncResult) * (IAsyncResult -> 'T) * (IAsyncResult -> unit)
  static member AwaitEvent : event:IEvent<'Del,'T> * ?cancelAction:(unit -> unit) -> Async<'T> (requires delegate and 'Del :> Delegate)
  static member AwaitIAsyncResult : iar:IAsyncResult * ?millisecondsTimeout:int -> Async<bool>
  static member AwaitTask : task:Task -> Async<unit>
  static member AwaitTask : task:Task<'T> -> Async<'T>
  static member AwaitWaitHandle : waitHandle:WaitHandle * ?millisecondsTimeout:int -> Async<bool>
  static member CancelDefaultToken : unit -> unit
  static member Catch : computation:Async<'T> -> Async<Choice<'T,exn>>
  static member Choice : computations:seq<Async<'T option>> -> Async<'T option>
  static member FromBeginEnd : beginAction:(AsyncCallback * obj -> IAsyncResult) * endAction:(IAsyncResult -> 'T) * ?cancelAction:(unit -> unit) -> Async<'T>
  ...

--------------------
type Async<'T> =
static member Async.RunSynchronously : computation:Async<'T> * ?timeout:int * ?cancellationToken:System.Threading.CancellationToken -> 'T
val rat : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo
val human : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo
val sheep : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo
val proteinInfos : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo []
val sequences : BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> []
Multiple items
module Array

from FSharpAux

--------------------
module Array

from Microsoft.FSharp.Collections
val map : mapping:('T -> 'U) -> array:'T [] -> 'U []
val x : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo
property EbiAPI.ProteinsAPIschema.ProteinFeatureInfo.Sequence: string with get, set
module BioArray

from BioFSharp
val ofAminoAcidSymbolString : s:#seq<char> -> BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol>
val names : string []
property EbiAPI.ProteinsAPIschema.ProteinFeatureInfo.EntryName: string with get, set
val align : a:BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> -> b:BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> -> Alignment.Alignment<AminoAcidSymbols.AminoAcidSymbol list,PairwiseAlignment.Score>


Creates the pairwise global alignment of two amino acid sequences based on the BLOSUM62 scoring matrix
val a : BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol>
val b : BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol>
val scoring : (AminoAcidSymbols.AminoAcidSymbol -> AminoAcidSymbols.AminoAcidSymbol -> int)
module ScoringMatrix

from BioFSharp.Algorithm
val getScoringMatrixAminoAcid : scoringMatrixType:ScoringMatrix.ScoringMatrixAminoAcid -> (AminoAcidSymbols.AminoAcidSymbol -> AminoAcidSymbols.AminoAcidSymbol -> int)
union case ScoringMatrix.ScoringMatrixAminoAcid.BLOSUM62: ScoringMatrix.ScoringMatrixAminoAcid
val costs : PairwiseAlignment.Costs<AminoAcidSymbols.AminoAcidSymbol>
module PairwiseAlignment

from BioFSharp.Algorithm
type Costs<'a> =
  { Open: int
    Continuation: int
    Similarity: 'a -> 'a -> int }
module AminoAcidSymbols

from BioFSharp
type AminoAcidSymbol =
  struct
    interface IBioItem
    interface IComparable
    private new : value:byte -> AminoAcidSymbol
    override Equals : other:obj -> bool
    override GetHashCode : unit -> int
    override ToString : unit -> string
    member private Value : byte
    static member Ala : AminoAcidSymbol
    static member Arg : AminoAcidSymbol
    static member Asn : AminoAcidSymbol
    ...
  end
module NeedlemanWunsch

from BioFSharp.Algorithm.PairwiseAlignment
val runAminoAcidSymbol : costs:PairwiseAlignment.Costs<AminoAcidSymbols.AminoAcidSymbol> -> fstSeq:BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> -> sndSeq:BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> -> Alignment.Alignment<AminoAcidSymbols.AminoAcidSymbol list,PairwiseAlignment.Score>
val createAlignmentMatrix : sequences:BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> [] -> float [,]


Aligns two amino acid sequences and writes their relative alignment score into a square matrix
type BioArray<'a (requires 'a :> IBioItem)> = 'a array
val mat : float [,]
Multiple items
module Array2D

from FSharpAux

--------------------
module Array2D

from Microsoft.FSharp.Collections
val zeroCreate : length1:int -> length2:int -> 'T [,]
val i : int
val j : int
val alignment : Alignment.Alignment<AminoAcidSymbols.AminoAcidSymbol list,PairwiseAlignment.Score>
Multiple items
val float : value:'T -> float (requires member op_Explicit)

--------------------
type float = System.Double

--------------------
type float<'Measure> = float
Alignment.Alignment.MetaData: PairwiseAlignment.Score
val mat' : float [,]
val map : mapping:('T -> 'U) -> array:'T [,] -> 'U [,]
val x : float
val maxBy : projection:('a -> 'b) -> arr:'a [,] -> 'a (requires comparison)
val id : x:'T -> 'T
val heatMap : GenericChart.GenericChart
val toJaggedArray : arr:'a [,] -> 'a [] []
val x : float [] []
type Chart =
  static member Area : xy:seq<#IConvertible * #IConvertible> * ?Name:string * ?ShowMarkers:bool * ?Showlegend:bool * ?MarkerSymbol:Symbol * ?Color:'a2 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font * ?Dash:DrawingStyle * ?Width:'a4 -> GenericChart
  static member Area : x:seq<#IConvertible> * y:seq<#IConvertible> * ?Name:string * ?ShowMarkers:bool * ?Showlegend:bool * ?MarkerSymbol:Symbol * ?Color:'a2 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font * ?Dash:DrawingStyle * ?Width:'a4 -> GenericChart
  static member Bar : keysvalues:seq<#IConvertible * #IConvertible> * ?Name:string * ?Showlegend:bool * ?Color:'a2 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font * ?Marker:Marker -> GenericChart
  static member Bar : keys:seq<#IConvertible> * values:seq<#IConvertible> * ?Name:string * ?Showlegend:bool * ?Color:'a2 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font * ?Marker:Marker -> GenericChart
  static member BoxPlot : xy:seq<'a0 * 'a1> * ?Name:string * ?Showlegend:bool * ?Color:'a2 * ?Fillcolor:'a3 * ?Opacity:float * ?Whiskerwidth:'a4 * ?Boxpoints:Boxpoints * ?Boxmean:BoxMean * ?Jitter:'a5 * ?Pointpos:'a6 * ?Orientation:Orientation -> GenericChart
  static member BoxPlot : ?x:'a0 * ?y:'a1 * ?Name:string * ?Showlegend:bool * ?Color:'a2 * ?Fillcolor:'a3 * ?Opacity:float * ?Whiskerwidth:'a4 * ?Boxpoints:Boxpoints * ?Boxmean:BoxMean * ?Jitter:'a5 * ?Pointpos:'a6 * ?Orientation:Orientation -> GenericChart
  static member Bubble : xysizes:seq<#IConvertible * #IConvertible * #IConvertible> * ?Name:string * ?Showlegend:bool * ?MarkerSymbol:Symbol * ?Color:'a3 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font -> GenericChart
  static member Bubble : x:seq<#IConvertible> * y:seq<#IConvertible> * sizes:seq<#IConvertible> * ?Name:string * ?Showlegend:bool * ?MarkerSymbol:Symbol * ?Color:'a3 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font -> GenericChart
  static member ChoroplethMap : locations:seq<string> * z:seq<#IConvertible> * ?Text:seq<#IConvertible> * ?Locationmode:LocationFormat * ?Autocolorscale:bool * ?Colorscale:Colorscale * ?Colorbar:'a2 * ?Marker:Marker * ?Zmin:'a3 * ?Zmax:'a4 -> GenericChart
  static member Column : keysvalues:seq<#IConvertible * #IConvertible> * ?Name:string * ?Showlegend:bool * ?Color:'a2 * ?Opacity:float * ?Labels:seq<#IConvertible> * ?TextPosition:TextPosition * ?TextFont:Font * ?Marker:Marker -> GenericChart
  ...
static member Chart.Heatmap : data:seq<#seq<'b>> * ?ColNames:seq<#System.IConvertible> * ?RowNames:seq<#System.IConvertible> * ?Name:string * ?Showlegend:bool * ?Opacity:float * ?Colorscale:StyleParam.Colorscale * ?Showscale:'e * ?Xgap:'f * ?Ygap:'g * ?zSmooth:StyleParam.SmoothAlg * ?Colorbar:'h -> GenericChart.GenericChart (requires 'b :> System.IConvertible)
module StyleParam

from FSharp.Plotly
type Colorscale =
  | Custom of seq<float * string>
  | RdBu
  | Earth
  | Blackbody
  | YIOrRd
  | YIGnBu
  | Bluered
  | Portland
  | Electric
  | Jet
  ...
    static member convert : (Colorscale -> obj)
union case StyleParam.Colorscale.Custom: seq<float * string> -> StyleParam.Colorscale
type Location =
  | Intracellular
  | Transmembrane
  | Extracellular
union case Location.Intracellular: Location
union case Location.Transmembrane: Location
union case Location.Extracellular: Location
val getLocationInfo : proteinInfo:EbiAPI.ProteinsAPIschema.ProteinFeatureInfo -> Location []


Maps protein information to information of location of amino acids of this protein
val proteinInfo : EbiAPI.ProteinsAPIschema.ProteinFeatureInfo
type ProteinsAPIschema =
  nested type Absorption
  nested type AlternativeName
  nested type AssociationType
  nested type CanonicalGeneType
  nested type CitationType
  nested type Client
  nested type CofactorType
  nested type CommentType
  nested type Component
  nested type ComponentType
  nested type Conflict
  nested type DbReference
  nested type DbReferenceObject
  nested type DbReferenceType
  nested type Disease
  nested type DiseaseComment
  nested type Domain
  nested type Entry
  nested type ErrorMessage
  nested type EventType
  nested type Evidence
  nested type EvidenceType
  nested type EvidencedString
  nested type EvidencedStringType
  nested type ExonMapType
  nested type Feature
  nested type FeatureType
  nested type GeneLocationType
  nested type GeneNameType
  nested type GeneType
  nested type GenomicLocationType
  nested type GnCoordinateType
  nested type GnEntry
  nested type GnFeatureType
  nested type ImportedFromType
  nested type IntAct
  nested type InteractantType
  nested type IsoformType
  nested type JournalType
  nested type KeywordType
  nested type Kinetics
  nested type Lineage
  nested type Link
  nested type LocationType
  nested type MoleculeType
  nested type Name
  nested type NameListType
  nested type OrganismNameType
  nested type OrganismType
  nested type PLocation2GLocation
  nested type PLocation2GLocationCollection
  nested type PeptideType
  nested type PhDependence
  nested type PolyphenPredictionType
  nested type PositionType
  nested type PropertyType
  nested type ProteinExistenceType
  nested type ProteinFeatureInfo
  nested type ProteinNameType
  nested type ProteinType
  nested type Proteome
  nested type RecommendedName
  nested type RedoxPotential
  nested type RedundantProteomeType
  nested type ReferenceType
  nested type Sequence
  nested type SequenceType1
  nested type SiftPredictionType
  nested type SourceDataType
  nested type SourceType
  nested type StatusType
  nested type SubcellLocationComment
  nested type SubcellularLocation
  nested type SubcellularLocationType
  nested type SubmissionType
  nested type SubmittedName
  nested type TemperatureDependence
  nested type UPInteraction
  nested type VariantType
type ProteinFeatureInfo =
  new : unit -> ProteinFeatureInfo + 1 overload
  member Accession : string with get, set
  member EntryName : string with get, set
  member Features : Feature[] with get, set
  member GeteGeneId : string with get, set
  member GeteProteinId : string with get, set
  member Sequence : string with get, set
  member SequenceChecksum : string with get, set
  member Taxid : FSharpOption<int> with get, set
  member ToString : unit -> string
  ...
property EbiAPI.ProteinsAPIschema.ProteinFeatureInfo.Features: EbiAPI.ProteinsAPIschema.Feature [] with get, set
val choose : chooser:('T -> 'U option) -> array:'T [] -> 'U []
val x : EbiAPI.ProteinsAPIschema.Feature
property EbiAPI.ProteinsAPIschema.Feature.Type: string with get, set
union case Option.Some: Value: 'T -> Option<'T>
property EbiAPI.ProteinsAPIschema.Feature.Begin: string with get, set
property EbiAPI.ProteinsAPIschema.Feature.End: string with get, set
property EbiAPI.ProteinsAPIschema.Feature.Description: string with get, set
union case Option.None: Option<'T>
val collect : mapping:('T -> 'U []) -> array:'T [] -> 'U []
val s : string
val e : string
val info : string
val init : count:int -> initializer:(int -> 'T) -> 'T []
Multiple items
val int : value:'T -> int (requires member op_Explicit)

--------------------
type int = int32

--------------------
type int<'Measure> = int
val getHydrophobicity : a:BioArray.BioArray<AminoAcidSymbols.AminoAcidSymbol> -> float []


Evaluates the hydrophobicity of all amino acids of a sequence
val hydro : (AminoAcidSymbols.AminoAcidSymbol -> float)
module AminoProperties

from BioFSharp
val initGetAminoProperty : property:AminoProperties.AminoProperty -> (AminoAcidSymbols.AminoAcidSymbol -> float)
union case AminoProperties.AminoProperty.HydrophobicityIndex: AminoProperties.AminoProperty
val ofWindowedBioArray : n:int -> pf:('a -> float) -> source:BioArray.BioArray<'a> -> float [] (requires 'a :> IBioItem)
val createBackgroundShape : i:int * location:Location -> Shape option


Creates the squares depicting the cell membrane
val location : Location
val i : float
Multiple items
type Shape =
  inherit DynamicObj
  new : unit -> Shape
  static member init : ?ShapeType:ShapeType * ?X0:'a * ?X1:'b * ?Y0:'c * ?Y1:'d * ?Path:'e * ?Opacity:'f * ?Line:Line * ?Fillcolor:'g * ?Layer:Layer * ?Xref:'h * ?Yref:'i -> Shape
  static member style : ?ShapeType:ShapeType * ?X0:'a0 * ?X1:'a1 * ?Y0:'a2 * ?Y1:'a3 * ?Path:'a4 * ?Opacity:'a5 * ?Line:Line * ?Fillcolor:'a6 * ?Layer:Layer * ?Xref:'a7 * ?Yref:'a8 -> (Shape -> Shape)

--------------------
new : unit -> Shape
static member Shape.init : ?ShapeType:StyleParam.ShapeType * ?X0:'a * ?X1:'b * ?Y0:'c * ?Y1:'d * ?Path:'e * ?Opacity:'f * ?Line:Line * ?Fillcolor:'g * ?Layer:StyleParam.Layer * ?Xref:'h * ?Yref:'i -> Shape
type ShapeType =
  | Circle
  | Rectangle
  | SvgPath
  | Line
    static member convert : (ShapeType -> obj)
    static member toString : (ShapeType -> string)
union case StyleParam.ShapeType.Rectangle: StyleParam.ShapeType
Multiple items
module Colors

from FSharpAux

--------------------
module Colors

from FSharp.Plotly
Multiple items
val toWebColor : c:Color -> string

--------------------
val toWebColor : c:Colors.Color -> string
Multiple items
module Table

from FSharpAux.Colors

--------------------
module Table

from FSharp.Plotly.Colors
Multiple items
module Office

from FSharpAux.Colors.Table

--------------------
module Office

from FSharp.Plotly.Colors.Table
Multiple items
val orange : Color

--------------------
val orange : Colors.Color
val hydrophobicityArr : float []
val mouseLocInfo : Location []
val shapes : GenericChart.GenericChart
static member Chart.Line : xy:seq<#System.IConvertible * #System.IConvertible> * ?Name:string * ?ShowMarkers:bool * ?Showlegend:bool * ?MarkerSymbol:StyleParam.Symbol * ?Color:'c * ?Opacity:float * ?Labels:seq<#System.IConvertible> * ?TextPosition:StyleParam.TextPosition * ?TextFont:Font * ?Dash:'e * ?Width:'f -> GenericChart.GenericChart
static member Chart.Line : x:seq<#System.IConvertible> * y:seq<#System.IConvertible> * ?Name:string * ?ShowMarkers:bool * ?Showlegend:bool * ?MarkerSymbol:StyleParam.Symbol * ?Color:'a2 * ?Opacity:float * ?Labels:seq<#System.IConvertible> * ?TextPosition:StyleParam.TextPosition * ?TextFont:Font * ?Dash:'a4 * ?Width:'a5 -> GenericChart.GenericChart
val indexed : array:'T [] -> (int * 'T) []
type Color =
  { A: byte
    R: byte
    G: byte
    B: byte }
Multiple items
val toHex : prefix:bool -> c:Color -> string

--------------------
val toHex : prefix:bool -> c:Colors.Color -> string
Multiple items
val blue : Color

--------------------
val blue : Colors.Color
static member Chart.withShapes : shapes:seq<Shape> -> (GenericChart.GenericChart -> GenericChart.GenericChart)
Multiple items
module Seq

from FSharpAux

--------------------
module Seq

from FSharp.Plotly

--------------------
module Seq

from Microsoft.FSharp.Collections
val indexed : source:seq<'T> -> seq<int * 'T>
val choose : chooser:('T -> 'U option) -> source:seq<'T> -> seq<'U>
static member Chart.withX_AxisStyle : title:'a * ?MinMax:(float * float) * ?Showgrid:'b * ?Showline:'c * ?Side:StyleParam.Side * ?Overlaying:StyleParam.AxisAnchorId * ?Id:int * ?Domain:(float * float) * ?Position:float * ?Anchor:StyleParam.AxisAnchorId -> (GenericChart.GenericChart -> GenericChart.GenericChart)
static member Chart.withY_AxisStyle : title:'a * ?MinMax:(float * float) * ?Showgrid:'b * ?Showline:'c * ?Side:StyleParam.Side * ?Overlaying:StyleParam.AxisAnchorId * ?Id:int * ?Domain:(float * float) * ?Position:float * ?Anchor:StyleParam.AxisAnchorId -> (GenericChart.GenericChart -> GenericChart.GenericChart)
Fork me on GitHub