ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_container::assignment(=) Interface Reference

Assign the contents of the input source to the output destin. More...

Detailed Description

Assign the contents of the input source to the output destin.

Parameters
[out]destin: The output scalar or array of the arbitrary rank and shape that can be of,
  1. type css_type
  2. type csi_type
  3. type csl_type
  4. type csc_type
  5. type csr_type
  6. type css_pdt
  7. type csi_pdt
  8. type csl_pdt
  9. type csc_pdt
  10. type csr_pdt
whose contents will be the taken from the contents of the input source argument upon return.
[in]source: The input scalar, or array of the same rank as the input array-like destin, of the same type and kind as destin, whose contents will be assigned to the output destin.


Possible calling interfaces

use pm_kind, only: LK
use pm_container, only: assignment(=)
destin = source
destin(..) = source
destin(..) = source(..)
This module contains the derived types for generating allocatable containers of scalar,...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
Warning
The condition rank(source) == 0 .or. rank(destin) == rank(source) must hold (i.e., only elemental assignments are possible).
Remarks
The procedures under discussion are pure.
The procedures under discussion are elemental.
This generic primarily exists to circumvent the gfortran bug as of version 11 in intrinsic assignment of arrays of type css_pdt.
See also
operator(==)
operator(<)
operator(>)
operator(>=)
operator(<=)
operator(/=)
assignment(=)


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_io, only: display_type
5 use pm_container, only: assignment(=)
6
7 implicit none
8
9 block
10 use pm_kind, only: SKG => SK
11 use pm_container, only: css_type
12 type(css_type) :: destin(3)
13
14 type(display_type) :: disp
15 disp = display_type(file = "main.out.F90")
16
17 call disp%skip()
18 call disp%show("destin(1) = css_type(SKG_'ParaMonte')")
19 destin(1) = css_type(SKG_'ParaMonte')
20 call disp%show("destin(1)")
21 call disp%show( destin(1) , deliml = SK_"""" )
22 call disp%skip()
23
24 call disp%skip()
25 call disp%show("destin = css_type(SKG_'ParaMonte')")
26 destin = css_type(SKG_'ParaMonte')
27 call disp%show("destin")
28 call disp%show( destin , deliml = SK_"""" )
29 call disp%skip()
30
31 call disp%skip()
32 call disp%show("destin = [css_type(SKG_'a'), css_type(SKG_'b'), css_type(SKG_'c')]")
33 destin = [css_type(SKG_'a'), css_type(SKG_'b'), css_type(SKG_'c')]
34 call disp%show("destin")
35 call disp%show( destin , deliml = SK_"""" )
36 call disp%skip()
37 end block
38
39#if PDT_ENABLED
40 block
41 use pm_kind, only: SKG => SK
42 use pm_container, only: css_pdt
43 type(css_pdt(SKG)) :: destin(3)
44
45 type(display_type) :: disp
46 disp = display_type(file = "main.out.F90")
47
48 call disp%skip()
49 call disp%show("destin(1) = css_pdt(SKG_'ParaMonte')")
50 destin(1) = css_pdt(SKG_'ParaMonte')
51 call disp%show("destin(1)")
52 call disp%show( destin(1) , deliml = SK_"""" )
53 call disp%skip()
54
55 call disp%skip()
56 call disp%show("destin = css_pdt(SKG_'ParaMonte')")
57 destin = css_pdt(SKG_'ParaMonte')
58 call disp%show("destin")
59 call disp%show( destin , deliml = SK_"""" )
60 call disp%skip()
61
62 call disp%skip()
63 call disp%show("destin = [css_pdt(SKG_'a'), css_pdt(SKG_'b'), css_pdt(SKG_'c')]")
64 destin = [css_pdt(SKG_'a'), css_pdt(SKG_'b'), css_pdt(SKG_'c')]
65 call disp%show("destin")
66 call disp%show( destin , deliml = SK_"""" )
67 call disp%skip()
68 end block
69#endif
70
71end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
This is the css_pdt parameterized type for generating instances of container of scalar of string obje...
This is the css_type type for generating instances of container of scalar of string objects.
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2destin(1) = css_type(SKG_'ParaMonte')
3destin(1)
4"ParaMonte"
5
6
7destin = css_type(SKG_'ParaMonte')
8destin
9"ParaMonte", "ParaMonte", "ParaMonte"
10
11
12destin = [css_type(SKG_'a'), css_type(SKG_'b'), css_type(SKG_'c')]
13destin
14"a", "b", "c"
15
16
Test:
test_pm_container
Bug:

Status: Unresolved
Source: GNU Fortran Compiler gfortran version 10.3-11
Description: The elemental implementations of the procedures under this generic interface yield incorrect results with gfortran 10.3.
Remedy (as of ParaMonte Library version 2.0.0): Currently unknown.
Todo:
Very Low Priority: The functionality of this generic interface can be extended to input arrays of higher rank.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, April 21, 2017, 3:54 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin

Definition at line 4227 of file pm_container.F90.


The documentation for this interface was generated from the following file: