Here is the code that does what I explained. I just switched up step 5 with 
a much simpler approach. Note that you can simply make this into a 
function/class template so that you could read general tensors, not just 
3x3 tensors. I leave that to you :-)

#include <deal.II/base/tensor.h>
#include <deal.II/base/parameter_handler.h>
#include <deal.II/base/utilities.h>
#include <string>
#include <iostream>
#include <sstream>

using namespace dealii;

int main() {


    // Let's read this tensor from an istringstream. You may read it from 
istream in your
    // application
    std::istringstream is("set Tensor3x3 = 11, 12, 13, 21, 22, 23, 31, 32, 33");


    // Now let's set up the ParameterHandler object. Since you are reading 3x3 
tensor, you are
    // dealing with 9 components and the default value should also have 9 
components (here 9 
    // zeros) in it
    ParameterHandler prm;

    using Patterns::List;
    using Patterns::Double;
    prm.declare_entry("Tensor3x3",
                      "0,0,0,0,0,0,0,0,0",
                      List(Double(), 9, 9),
                      "Tensor components");


    // Let's read it in now
    prm.parse_input(is);
    const std::string tensorString = prm.get("Tensor3x3");

    // Now let's split the string into individual components
    const std::vector<std::string> tensorComponentsString =
            Utilities::split_string_list(tensorString);

    // Now go over each component, convert it into double and put it in the 
tensor
    Tensor<2, 3> tensor;
    unsigned int counter = 0;
    for (unsigned ii = 0; ii < tensor.dimension; ++ii) {
        for (unsigned jj = 0; jj < tensor.dimension; ++jj) {
            tensor[ii][jj] = 
Utilities::string_to_double(tensorComponentsString[counter++]);
        }
    }


    // Let's check if we read the tensor correctly
    for (unsigned ii = 0; ii < tensor.dimension; ++ii) {
        for (unsigned jj = 0; jj < tensor.dimension; ++jj) {
            std::cout << "component [" << ii << ", " << jj << " ] = "
                      << tensor[ii][jj] << std::endl;
        }
    }


}






On Wednesday, April 15, 2020 at 3:54:57 AM UTC-4, Paras Kumar wrote:
>
> Ahmad,
>
> On Tuesday, April 14, 2020 at 5:39:06 PM UTC+2, Ahmad Shahba wrote:
>>
>> I don't know if it is the optimal way but I would use the following 
>> approach
>>
>>    1. Read tensor components as Patterns::List 
>>    2. Use  get 
>>    
>> <https://www.dealii.org/current/doxygen/deal.II/classParameterHandler.html#a91cfbaca954f444047302446a4e87125>
>>   
>>    method of ParameterHandler to read all tensor  components as one string 
>>    into a string variable
>>    3. Use  split_string_list 
>>    
>> <https://www.dealii.org/current/doxygen/deal.II/namespaceUtilities.html#a0a2b9ceded96f25b47eccd7c177dc67e>
>>       in the Utilities namespace to split the string into individual 
>>    components (still all components are strings)
>>    4. Convert  string-type components into doubles using  
>>    string_to_double 
>>    
>> <https://www.dealii.org/current/doxygen/deal.II/namespaceUtilities.html#ab3177021843ad87857e6f8c5d98d29a5>
>>  in 
>>    the Utilities namespace
>>    5. Use  Tensor 
>>    
>> <https://www.dealii.org/current/doxygen/deal.II/classTensor.html#a6118e90e58ca4041ba452ab647db0b29>
>>     (const Tensor 
>>    <https://www.dealii.org/current/doxygen/deal.II/classTensor.html>< 
>>    rank_, dim, OtherNumber > &initializer) to construct the tensor
>>
>> Could you please elaborate on pointa 1 and 5, or share a code snippet 
> which explains the process. 
>
> Best,
> Paras
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/18dbabf5-84a8-48a6-885e-b955f0566350%40googlegroups.com.

Reply via email to