DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_STKMeshField_impl.hpp
1 //---------------------------------------------------------------------------//
2 /*
3  Copyright (c) 2012, Stuart R. Slattery
4  All rights reserved.
5 
6  Redistribution and use in source and binary forms, with or without
7  modification, are permitted provided that the following conditions are
8  met:
9 
10  *: Redistributions of source code must retain the above copyright
11  notice, this list of conditions and the following disclaimer.
12 
13  *: Redistributions in binary form must reproduce the above copyright
14  notice, this list of conditions and the following disclaimer in the
15  documentation and/or other materials provided with the distribution.
16 
17  *: Neither the name of the University of Wisconsin - Madison nor the
18  names of its contributors may be used to endorse or promote products
19  derived from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
24  A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25  HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26  SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
27  LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
28  DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
29  THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
31  OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 */
33 //---------------------------------------------------------------------------//
39 //---------------------------------------------------------------------------//
40 
41 #ifndef DTK_STKMESHFIELD_IMPL_HPP
42 #define DTK_STKMESHFIELD_IMPL_HPP
43 
44 #include <vector>
45 
46 #include "DTK_DBC.hpp"
47 
48 #include <Teuchos_DefaultMpiComm.hpp>
49 
50 #include <stk_mesh/base/FieldParallel.hpp>
51 #include <stk_mesh/base/FieldRestriction.hpp>
52 #include <stk_mesh/base/GetEntities.hpp>
53 #include <stk_mesh/base/MetaData.hpp>
54 #include <stk_mesh/base/Selector.hpp>
55 
56 namespace DataTransferKit
57 {
58 //---------------------------------------------------------------------------//
59 // Constructor.
60 template <class Scalar, class FieldType>
62  const Teuchos::RCP<stk::mesh::BulkData> &bulk_data,
63  const Teuchos::Ptr<FieldType> &field, const int field_dim )
64  : d_bulk_data( bulk_data )
65  , d_field( field )
66  , d_field_dim( field_dim )
67 {
68  // Get the field restriction.
69  const stk::mesh::FieldRestrictionVector field_restrictions =
70  d_field->restrictions();
71 
72  // Build a composite selector.
73  stk::mesh::Selector field_selector;
74  for ( stk::mesh::FieldRestriction r : field_restrictions )
75  {
76  field_selector = field_selector | r.selector();
77  }
78 
79  // Get the buckets for the field entity rank.
80  const stk::mesh::BucketVector &field_buckets =
81  field_selector.get_buckets( d_field->entity_rank() );
82 
83  // Get the local entities over which the field is defined.
84  stk::mesh::get_selected_entities( field_selector, field_buckets,
85  d_field_entities );
86 
87  // Map the field entity ids.
88  int num_entities = d_field_entities.size();
89  for ( int n = 0; n < num_entities; ++n )
90  {
91  d_id_map.emplace( d_bulk_data->identifier( d_field_entities[n] ), n );
92  }
93 
94  // Get the locally owned field entity ids.
95  for ( int n = 0; n < num_entities; ++n )
96  {
97  if ( d_bulk_data->parallel_owner_rank( d_field_entities[n] ) ==
98  d_bulk_data->parallel_rank() )
99  {
100  d_support_ids.push_back(
101  d_bulk_data->identifier( d_field_entities[n] ) );
102  }
103  }
104 }
105 
106 //---------------------------------------------------------------------------//
107 // Get the dimension of the field.
108 template <class Scalar, class FieldType>
110 {
111  return d_field_dim;
112 }
113 
114 //---------------------------------------------------------------------------//
115 // Get the locally-owned entity support location ids of the field.
116 template <class Scalar, class FieldType>
117 Teuchos::ArrayView<const SupportId>
119 {
120  return d_support_ids();
121 }
122 
123 //---------------------------------------------------------------------------//
124 // Given a local support id and a dimension, read data from the application
125 // field.
126 template <class Scalar, class FieldType>
127 double
129  const int dimension ) const
130 {
131  DTK_REQUIRE( d_id_map.count( support_id ) );
132  int local_id = d_id_map.find( support_id )->second;
133  return stk::mesh::field_data( *d_field,
134  d_field_entities[local_id] )[dimension];
135 }
136 
137 //---------------------------------------------------------------------------//
138 // Given a local support id, dimension, and field value, write data into the
139 // application field.
140 template <class Scalar, class FieldType>
142  const SupportId support_id, const int dimension, const double data )
143 {
144  DTK_REQUIRE( d_id_map.count( support_id ) );
145  int local_id = d_id_map.find( support_id )->second;
146  stk::mesh::field_data( *d_field, d_field_entities[local_id] )[dimension] =
147  data;
148 }
149 
150 //---------------------------------------------------------------------------//
151 // Finalize a field after writing into it.
152 template <class Scalar, class FieldType>
154 {
155  stk::mesh::copy_owned_to_shared(
156  *d_bulk_data,
157  std::vector<const stk::mesh::FieldBase *>( 1, d_field.getRawPtr() ) );
158 }
159 
160 //---------------------------------------------------------------------------//
161 
162 } // end namespace DataTransferKit
163 
164 //---------------------------------------------------------------------------//
165 
166 #endif // end DTK_STKMESHFIELD_IMPL_HPP
167 
168 //---------------------------------------------------------------------------//
169 // end DTK_STKMeshField_impl.hpp
170 //---------------------------------------------------------------------------//
double readFieldData(const SupportId support_id, const int dimension) const override
Given a local support id and a dimension, read data from the application field.
STKMeshField(const Teuchos::RCP< stk::mesh::BulkData > &bulk_data, const Teuchos::Ptr< FieldType > &field, const int field_dim)
Constructor.
unsigned long int SupportId
Support id type.
Definition: DTK_Types.hpp:57
void finalizeAfterWrite() override
Finalize a field after writing into it.
Assertions and Design-by-Contract for error handling.
void writeFieldData(const SupportId support_id, const int dimension, const double data) override
Given a local support id, dimension, and field value, write data into the application field...
int dimension() const override
Get the dimension of the field.
DTK_BasicEntitySet.cpp.
Teuchos::ArrayView< const SupportId > getLocalSupportIds() const override
Get the locally-owned entity support location ids of the field.